Belle II Software  release-08-00-10
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 <TDirectory.h>
18 #include <TFile.h>
19 #include <TH1F.h>
20 #include <TMath.h>
21 #include <TString.h>
22 
23 /* C++ headers. */
24 #include <iostream>
25 
26 using namespace std;
27 using namespace Belle2;
28 using namespace ECL;
29 using namespace Calibration;
30 
31 
33 eclEdgeAlgorithm::eclEdgeAlgorithm(): CalibrationAlgorithm("eclEdgeCollector")
34 {
36  "Generate payloads ECLCrystalThetaEdge and ECLCrystalPhiEdge found by eclEdgeCollector"
37  );
38 }
39 
40 
43 {
44 
45  //-----------------------------------------------------------------------------------
46  //..Read in histograms created by the collector
47  auto eclCrystalX = getObjectPtr<TH1F>("eclCrystalX");
48  auto eclCrystalY = getObjectPtr<TH1F>("eclCrystalY");
49  auto eclCrystalZ = getObjectPtr<TH1F>("eclCrystalZ");
50  auto eclCrystalR = getObjectPtr<TH1F>("eclCrystalR");
51  auto eclCrystalTheta = getObjectPtr<TH1F>("eclCrystalTheta");
52  auto eclCrystalPhi = getObjectPtr<TH1F>("eclCrystalPhi");
53  auto eclCrystalDirTheta = getObjectPtr<TH1F>("eclCrystalDirTheta");
54  auto eclCrystalDirPhi = getObjectPtr<TH1F>("eclCrystalDirPhi");
55  auto eclCrystalEdgeTheta = getObjectPtr<TH1F>("eclCrystalEdgeTheta");
56  auto eclCrystalEdgePhi = getObjectPtr<TH1F>("eclCrystalEdgePhi");
57  auto eclEdgeCounter = getObjectPtr<TH1F>("eclEdgeCounter");
58 
59  //..And write them to disk
60  TFile* histfile = new TFile("eclEdgeAlgorithm.root", "recreate");
61  eclCrystalX->Write();
62  eclCrystalY->Write();
63  eclCrystalZ->Write();
64  eclCrystalR->Write();
65  eclCrystalTheta->Write();
66  eclCrystalPhi->Write();
67  eclCrystalDirTheta->Write();
68  eclCrystalDirPhi->Write();
69  eclCrystalEdgeTheta->Write();
70  eclCrystalEdgePhi->Write();
71  eclEdgeCounter->Write();
72 
73  //-----------------------------------------------------------------------------------
74  //..Set up
75 
76  //..Number of collector calls. Intended to be 1.
77  const int nCalls = (int)(eclEdgeCounter->GetBinContent(1) + 0.0001);
78 
79  //..First crysID of each thetaID will be handy
80  int firstCrysID[69] = {};
81  for (int thetaID = 1; thetaID < 69; thetaID++) {
82  firstCrysID[thetaID] = firstCrysID[thetaID - 1] + m_crystalsPerRing[thetaID - 1];
83  }
84 
85  //..Histograms of width payloads
86  TH1F* eclCrystalWidthTheta = new TH1F("eclCrystalWidthTheta", "Width of each crystal in theta;cellID;crystal width (rad)",
88  8737);
89  TH1F* eclCrystalWidthPhi = new TH1F("eclCrystalWidthPhi", "Width of each crystal in phi;cellID;crystal width (rad)",
91 
92  //-----------------------------------------------------------------------------------
93  //..Edges of the crystals in theta
94  std::vector<float> tempThetaEdge;
95  std::vector<float> tempUnc(ECLElementNumbers::c_NCrystals, 0.);
96  for (int cellID = 1; cellID <= ECLElementNumbers::c_NCrystals; cellID++) {
97  tempThetaEdge.push_back(eclCrystalEdgeTheta->GetBinContent(cellID) / nCalls);
98  }
99  ECLCrystalCalib* crystalThetaEdge = new ECLCrystalCalib();
100  crystalThetaEdge->setCalibVector(tempThetaEdge, tempUnc);
101  saveCalibration(crystalThetaEdge, "ECLCrystalThetaEdge");
102  B2RESULT("eclEdgeAlgorithm: successfully stored payload ECLCrystalThetaEdge");
103 
104  //-----------------------------------------------------------------------------------
105  //..Edges of the crystals in phi
106  std::vector<float> tempPhiEdge;
107  for (int cellID = 1; cellID <= ECLElementNumbers::c_NCrystals; cellID++) {
108  tempPhiEdge.push_back(eclCrystalEdgePhi->GetBinContent(cellID) / nCalls);
109  }
110  ECLCrystalCalib* crystalPhiEdge = new ECLCrystalCalib();
111  crystalPhiEdge->setCalibVector(tempPhiEdge, tempUnc);
112  saveCalibration(crystalPhiEdge, "ECLCrystalPhiEdge");
113  B2RESULT("eclEdgeAlgorithm: successfully stored payload ECLCrystalPhiEdge");
114 
115  //-----------------------------------------------------------------------------------
116  //..Width of the crystals in phi.
117  // Lower edge of next crystal in phi minus lower edge of crystal.
118  std::vector<float> tempPhiWidth;
119  int crysID = -1;
120  for (int thetaID = 0; thetaID < 69; thetaID++) {
121  for (int phiID = 0; phiID < m_crystalsPerRing[thetaID]; phiID++) {
122  crysID++;
123  int nextID = crysID + 1;
124  if (phiID == m_crystalsPerRing[thetaID] - 1) {nextID -= m_crystalsPerRing[thetaID];}
125  double width = tempPhiEdge[nextID] - tempPhiEdge[crysID];
126  if (width < 0) {width += 2.*TMath::Pi();}
127  tempPhiWidth.push_back(width);
128  }
129  }
130  ECLCrystalCalib* crystalPhiWidth = new ECLCrystalCalib();
131  crystalPhiWidth->setCalibVector(tempPhiWidth, tempUnc);
132  saveCalibration(crystalPhiWidth, "ECLCrystalPhiWidth");
133  B2RESULT("eclEdgeAlgorithm: successfully stored payload ECLCrystalPhiWidth");
134 
135  //..Also store in histogram
136  for (int cellID = 1; cellID <= ECLElementNumbers::c_NCrystals; cellID++) {
137  eclCrystalWidthPhi->SetBinContent(cellID, tempPhiWidth[cellID - 1]);
138  eclCrystalWidthPhi->SetBinError(cellID, 0.);
139  }
140  histfile->cd();
141  eclCrystalWidthPhi->Write();
142 
143  //-----------------------------------------------------------------------------------
144  //..Width in theta. Look crystals in the next thetaID that overlap in phi.
145  // Last thetaID is a special case.
146  std::vector<float> tempThetaWidth;
147  for (int thetaID = 0; thetaID < 68; thetaID++) {
148  for (int ic = firstCrysID[thetaID]; ic < firstCrysID[thetaID] + m_crystalsPerRing[thetaID]; ic++) {
149  double minThetaWidth = 999.;
150  double maxThetaWidth = -999.;
151  for (int icnext = firstCrysID[thetaID + 1]; icnext < firstCrysID[thetaID + 1] + m_crystalsPerRing[thetaID + 1]; icnext++) {
152 
153  //..Lower edge of ic falls within icnext
154  double offset = tempPhiEdge[ic] - tempPhiEdge[icnext];
155  if (offset < -TMath::Pi()) {offset += 2.*TMath::Pi();}
156  if (offset > TMath::Pi()) {offset -= 2.*TMath::Pi();}
157  if (offset >= 0. and offset < tempPhiWidth[icnext]) {
158  double width = tempThetaEdge[icnext] - tempThetaEdge[ic];
159  if (width < minThetaWidth) {minThetaWidth = width;}
160  if (width > maxThetaWidth) {maxThetaWidth = width;}
161  }
162 
163  //..Lower edge of icnext falls within ic
164  offset = tempPhiEdge[icnext] - tempPhiEdge[ic];
165  if (offset < -TMath::Pi()) {offset += 2.*TMath::Pi();}
166  if (offset > TMath::Pi()) {offset -= 2.*TMath::Pi();}
167  if (offset >= 0. and offset < tempPhiWidth[ic]) {
168  double width = tempThetaEdge[icnext] - tempThetaEdge[ic];
169  if (width < minThetaWidth) {minThetaWidth = width;}
170  if (width > maxThetaWidth) {maxThetaWidth = width;}
171  }
172 
173  }
174  tempThetaWidth.push_back(0.5 * (maxThetaWidth + minThetaWidth));
175  }
176  }
177 
178  //..Last thetaID; assume crystals end at nominal value from detector drawings
179  const double upperThetaEdge = 2.7416;
180  const int thetaID68 = 68;
181  for (int ic = firstCrysID[thetaID68]; ic < firstCrysID[thetaID68] + m_crystalsPerRing[thetaID68]; ic++) {
182  tempThetaWidth.push_back(upperThetaEdge - tempThetaEdge[ic]);
183  }
184 
185  //..Store the payload
186  ECLCrystalCalib* crystalThetaWidth = new ECLCrystalCalib();
187  crystalThetaWidth->setCalibVector(tempThetaWidth, tempUnc);
188  saveCalibration(crystalThetaWidth, "ECLCrystalThetaWidth");
189  B2RESULT("eclEdgeAlgorithm: successfully stored payload ECLCrystalThetaWidth");
190 
191  //..Also store in histogram
192  for (int cellID = 1; cellID <= ECLElementNumbers::c_NCrystals; cellID++) {
193  eclCrystalWidthTheta->SetBinContent(cellID, tempThetaWidth[cellID - 1]);
194  eclCrystalWidthTheta->SetBinError(cellID, 0.);
195  }
196  histfile->cd();
197  eclCrystalWidthTheta->Write();
198 
199  //-----------------------------------------------------------------------------------
200  //..Done
201  histfile->Close();
202  return c_OK;
203 }
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.