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