Belle II Software development
eclEdgeAlgorithm.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9/* Own header. */
10#include <ecl/calibration/eclEdgeAlgorithm.h>
11
12/* ECL headers. */
13#include <ecl/dataobjects/ECLElementNumbers.h>
14#include <ecl/dbobjects/ECLCrystalCalib.h>
15
16/* ROOT headers. */
17#include <TFile.h>
18#include <TH1F.h>
19#include <TMath.h>
20
21/* C++ headers. */
22#include <iostream>
23
24using namespace std;
25using namespace Belle2;
26using namespace ECL;
27using namespace Calibration;
28
29
32{
34 "Generate payloads ECLCrystalThetaEdge and ECLCrystalPhiEdge found by eclEdgeCollector"
35 );
36}
37
38
41{
42
43 //-----------------------------------------------------------------------------------
44 //..Read in histograms created by the collector
45 auto eclCrystalX = getObjectPtr<TH1F>("eclCrystalX");
46 auto eclCrystalY = getObjectPtr<TH1F>("eclCrystalY");
47 auto eclCrystalZ = getObjectPtr<TH1F>("eclCrystalZ");
48 auto eclCrystalR = getObjectPtr<TH1F>("eclCrystalR");
49 auto eclCrystalTheta = getObjectPtr<TH1F>("eclCrystalTheta");
50 auto eclCrystalPhi = getObjectPtr<TH1F>("eclCrystalPhi");
51 auto eclCrystalDirTheta = getObjectPtr<TH1F>("eclCrystalDirTheta");
52 auto eclCrystalDirPhi = getObjectPtr<TH1F>("eclCrystalDirPhi");
53 auto eclCrystalEdgeTheta = getObjectPtr<TH1F>("eclCrystalEdgeTheta");
54 auto eclCrystalEdgePhi = getObjectPtr<TH1F>("eclCrystalEdgePhi");
55 auto eclEdgeCounter = getObjectPtr<TH1F>("eclEdgeCounter");
56
57 //..And write them to disk
58 TFile* histfile = new TFile("eclEdgeAlgorithm.root", "recreate");
59 eclCrystalX->Write();
60 eclCrystalY->Write();
61 eclCrystalZ->Write();
62 eclCrystalR->Write();
63 eclCrystalTheta->Write();
64 eclCrystalPhi->Write();
65 eclCrystalDirTheta->Write();
66 eclCrystalDirPhi->Write();
67 eclCrystalEdgeTheta->Write();
68 eclCrystalEdgePhi->Write();
69 eclEdgeCounter->Write();
70
71 //-----------------------------------------------------------------------------------
72 //..Set up
73
74 //..Number of collector calls. Intended to be 1.
75 const int nCalls = (int)(eclEdgeCounter->GetBinContent(1) + 0.0001);
76
77 //..First crysID of each thetaID will be handy
78 int firstCrysID[69] = {};
79 for (int thetaID = 1; thetaID < 69; thetaID++) {
80 firstCrysID[thetaID] = firstCrysID[thetaID - 1] + m_crystalsPerRing[thetaID - 1];
81 }
82
83 //..Histograms of width payloads
84 TH1F* eclCrystalWidthTheta = new TH1F("eclCrystalWidthTheta", "Width of each crystal in theta;cellID;crystal width (rad)",
86 8737);
87 TH1F* eclCrystalWidthPhi = new TH1F("eclCrystalWidthPhi", "Width of each crystal in phi;cellID;crystal width (rad)",
89
90 //-----------------------------------------------------------------------------------
91 //..Edges of the crystals in theta
92 std::vector<float> tempThetaEdge;
93 std::vector<float> tempUnc(ECLElementNumbers::c_NCrystals, 0.);
94 for (int cellID = 1; cellID <= ECLElementNumbers::c_NCrystals; cellID++) {
95 tempThetaEdge.push_back(eclCrystalEdgeTheta->GetBinContent(cellID) / nCalls);
96 }
97 ECLCrystalCalib* crystalThetaEdge = new ECLCrystalCalib();
98 crystalThetaEdge->setCalibVector(tempThetaEdge, tempUnc);
99 saveCalibration(crystalThetaEdge, "ECLCrystalThetaEdge");
100 B2RESULT("eclEdgeAlgorithm: successfully stored payload ECLCrystalThetaEdge");
101
102 //-----------------------------------------------------------------------------------
103 //..Edges of the crystals in phi
104 std::vector<float> tempPhiEdge;
105 for (int cellID = 1; cellID <= ECLElementNumbers::c_NCrystals; cellID++) {
106 tempPhiEdge.push_back(eclCrystalEdgePhi->GetBinContent(cellID) / nCalls);
107 }
108 ECLCrystalCalib* crystalPhiEdge = new ECLCrystalCalib();
109 crystalPhiEdge->setCalibVector(tempPhiEdge, tempUnc);
110 saveCalibration(crystalPhiEdge, "ECLCrystalPhiEdge");
111 B2RESULT("eclEdgeAlgorithm: successfully stored payload ECLCrystalPhiEdge");
112
113 //-----------------------------------------------------------------------------------
114 //..Width of the crystals in phi.
115 // Lower edge of next crystal in phi minus lower edge of crystal.
116 std::vector<float> tempPhiWidth;
117 int crysID = -1;
118 for (int thetaID = 0; thetaID < 69; thetaID++) {
119 for (int phiID = 0; phiID < m_crystalsPerRing[thetaID]; phiID++) {
120 crysID++;
121 int nextID = crysID + 1;
122 if (phiID == m_crystalsPerRing[thetaID] - 1) {nextID -= m_crystalsPerRing[thetaID];}
123 double width = tempPhiEdge[nextID] - tempPhiEdge[crysID];
124 if (width < 0) {width += 2.*TMath::Pi();}
125 tempPhiWidth.push_back(width);
126 }
127 }
128 ECLCrystalCalib* crystalPhiWidth = new ECLCrystalCalib();
129 crystalPhiWidth->setCalibVector(tempPhiWidth, tempUnc);
130 saveCalibration(crystalPhiWidth, "ECLCrystalPhiWidth");
131 B2RESULT("eclEdgeAlgorithm: successfully stored payload ECLCrystalPhiWidth");
132
133 //..Also store in histogram
134 for (int cellID = 1; cellID <= ECLElementNumbers::c_NCrystals; cellID++) {
135 eclCrystalWidthPhi->SetBinContent(cellID, tempPhiWidth[cellID - 1]);
136 eclCrystalWidthPhi->SetBinError(cellID, 0.);
137 }
138 histfile->cd();
139 eclCrystalWidthPhi->Write();
140
141 //-----------------------------------------------------------------------------------
142 //..Width in theta. Look crystals in the next thetaID that overlap in phi.
143 // Last thetaID is a special case.
144 std::vector<float> tempThetaWidth;
145 for (int thetaID = 0; thetaID < 68; thetaID++) {
146 for (int ic = firstCrysID[thetaID]; ic < firstCrysID[thetaID] + m_crystalsPerRing[thetaID]; ic++) {
147 double minThetaWidth = 999.;
148 double maxThetaWidth = -999.;
149 for (int icnext = firstCrysID[thetaID + 1]; icnext < firstCrysID[thetaID + 1] + m_crystalsPerRing[thetaID + 1]; icnext++) {
150
151 //..Lower edge of ic falls within icnext
152 double offset = tempPhiEdge[ic] - tempPhiEdge[icnext];
153 if (offset < -TMath::Pi()) {offset += 2.*TMath::Pi();}
154 if (offset > TMath::Pi()) {offset -= 2.*TMath::Pi();}
155 if (offset >= 0. and offset < tempPhiWidth[icnext]) {
156 double width = tempThetaEdge[icnext] - tempThetaEdge[ic];
157 if (width < minThetaWidth) {minThetaWidth = width;}
158 if (width > maxThetaWidth) {maxThetaWidth = width;}
159 }
160
161 //..Lower edge of icnext falls within ic
162 offset = tempPhiEdge[icnext] - tempPhiEdge[ic];
163 if (offset < -TMath::Pi()) {offset += 2.*TMath::Pi();}
164 if (offset > TMath::Pi()) {offset -= 2.*TMath::Pi();}
165 if (offset >= 0. and offset < tempPhiWidth[ic]) {
166 double width = tempThetaEdge[icnext] - tempThetaEdge[ic];
167 if (width < minThetaWidth) {minThetaWidth = width;}
168 if (width > maxThetaWidth) {maxThetaWidth = width;}
169 }
170
171 }
172 tempThetaWidth.push_back(0.5 * (maxThetaWidth + minThetaWidth));
173 }
174 }
175
176 //..Last thetaID; assume crystals end at nominal value from detector drawings
177 const double upperThetaEdge = 2.7416;
178 const int thetaID68 = 68;
179 for (int ic = firstCrysID[thetaID68]; ic < firstCrysID[thetaID68] + m_crystalsPerRing[thetaID68]; ic++) {
180 tempThetaWidth.push_back(upperThetaEdge - tempThetaEdge[ic]);
181 }
182
183 //..Store the payload
184 ECLCrystalCalib* crystalThetaWidth = new ECLCrystalCalib();
185 crystalThetaWidth->setCalibVector(tempThetaWidth, tempUnc);
186 saveCalibration(crystalThetaWidth, "ECLCrystalThetaWidth");
187 B2RESULT("eclEdgeAlgorithm: successfully stored payload ECLCrystalThetaWidth");
188
189 //..Also store in histogram
190 for (int cellID = 1; cellID <= ECLElementNumbers::c_NCrystals; cellID++) {
191 eclCrystalWidthTheta->SetBinContent(cellID, tempThetaWidth[cellID - 1]);
192 eclCrystalWidthTheta->SetBinError(cellID, 0.);
193 }
194 histfile->cd();
195 eclCrystalWidthTheta->Write();
196
197 //-----------------------------------------------------------------------------------
198 //..Done
199 histfile->Close();
200 return c_OK;
201}
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
General DB object to store one calibration number per ECL crystal.
void setCalibVector(const std::vector< float > &CalibConst, const std::vector< float > &CalibConstUnc)
Set vector of constants with uncertainties.
const short m_crystalsPerRing[69]
crystals per thetaID
virtual EResult calibrate() override
..Run algorithm
const int c_NCrystals
Number of crystals.
Abstract base class for different kinds of events.
STL namespace.