Belle II Software  release-06-00-14
eclee5x5Algorithm.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 #include <ecl/calibration/eclee5x5Algorithm.h>
10 #include <ecl/dbobjects/ECLCrystalCalib.h>
11 
12 #include "TH2F.h"
13 #include "TFile.h"
14 #include "TF1.h"
15 #include "TROOT.h"
16 #include "TMatrixD.h"
17 #include "TMatrixDSym.h"
18 #include "TDecompLU.h"
19 
20 using namespace std;
21 using namespace Belle2;
22 using namespace ECL;
23 
25 eclee5x5Algorithm::eclee5x5Algorithm(): CalibrationAlgorithm("eclee5x5Collector")
26 {
28  "Perform energy calibration of ecl crystals by analyzing energy in 25-crystal sums from Bhabha events"
29  );
30 }
31 
33 {
36  TH1F* dummy;
37  dummy = (TH1F*)gROOT->FindObject("AverageExpECrys");
38  if (dummy) {delete dummy;}
39  dummy = (TH1F*)gROOT->FindObject("AverageElecCalib");
40  if (dummy) {delete dummy;}
41  dummy = (TH1F*)gROOT->FindObject("AverageInitCalib");
42  if (dummy) {delete dummy;}
43  dummy = (TH1F*)gROOT->FindObject("meanEnvsCrysID");
44  if (dummy) {delete dummy;}
45 
47  gROOT->SetBatch();
48 
51  B2INFO("eclee5x5Algorithm parameters:");
52  B2INFO("outputName = " << m_outputName);
53  B2INFO("minEntries = " << m_minEntries);
54  B2INFO("payloadName = " << m_payloadName);
55  B2INFO("storeConst = " << m_storeConst);
56  if (m_payloadName == "ECLeedPhiData" or m_payloadName == "ECLeedPhiMC" or m_payloadName == "None") {
57  B2INFO("fracLo = " << m_fracLo);
58  B2INFO("fracHiSym = " << m_fracHiSym);
59  B2INFO("fracHiASym = " << m_fracHiASym);
60  B2INFO("nsigLo = " << m_nsigLo);
61  B2INFO("nsigHiSym = " << m_nsigHiSym);
62  B2INFO("nsigHiASym = " << m_nsigHiASym);
63  }
64 
67  auto EnVsCrysID = getObjectPtr<TH2F>("EnVsCrysID");
68  auto RvsCrysID = getObjectPtr<TH1F>("RvsCrysID");
69  auto NRvsCrysID = getObjectPtr<TH1F>("NRvsCrysID");
70  auto Qmatrix = getObjectPtr<TH2F>("Qmatrix");
71  auto ElecCalibvsCrys = getObjectPtr<TH1F>("ElecCalibvsCrys");
72  auto ExpEvsCrys = getObjectPtr<TH1F>("ExpEvsCrys");
73  auto InitialCalibvsCrys = getObjectPtr<TH1F>("InitialCalibvsCrys");
74  auto CalibEntriesvsCrys = getObjectPtr<TH1F>("CalibEntriesvsCrys");
75  auto EntriesvsCrys = getObjectPtr<TH1F>("EntriesvsCrys");
76  auto dPhivsThetaID = getObjectPtr<TH2F>("dPhivsThetaID");
77 
81  TH1F* AverageExpECrys = new TH1F("AverageExpECrys", "Average expected E per crys from collector;Crystal ID;Energy (GeV)", 8736, 0,
82  8736);
83  TH1F* AverageElecCalib = new TH1F("AverageElecCalib", "Average electronics calib const vs crystal;Crystal ID;Calibration constant",
84  8736, 0, 8736);
85  TH1F* AverageInitCalib = new TH1F("AverageInitCalib", "Average initial calib const vs crystal;Crystal ID;Calibration constant",
86  8736, 0, 8736);
87  TH1F* meanEnvsCrysID = new TH1F("meanEnvsCrysID", "Mean normalized energy vs crystal;CrystalID;E/Eexp", 8736, 0, 8736);
88 
89  for (int cellID = 1; cellID <= 8736; cellID++) {
90  double TotEntries = CalibEntriesvsCrys->GetBinContent(cellID);
91  if (TotEntries > 0.) {
92  AverageElecCalib->SetBinContent(cellID, ElecCalibvsCrys->GetBinContent(cellID) / TotEntries);
93  AverageExpECrys->SetBinContent(cellID, ExpEvsCrys->GetBinContent(cellID) / TotEntries);
94  AverageInitCalib->SetBinContent(cellID, InitialCalibvsCrys->GetBinContent(cellID) / TotEntries);
95  }
96 
97  TH1D* En = EnVsCrysID->ProjectionY("En", cellID, cellID);
98  meanEnvsCrysID->SetBinContent(cellID, En->GetMean());
99  meanEnvsCrysID->SetBinError(cellID, En->GetStdDev());
100  }
101 
104  TString fName = m_outputName;
105  TFile* histfile = new TFile(fName, "recreate");
106  EnVsCrysID->Write();
107  RvsCrysID->Write();
108  NRvsCrysID->Write();
109  Qmatrix->Write();
110  AverageElecCalib->Write();
111  AverageExpECrys->Write();
112  AverageInitCalib->Write();
113  EntriesvsCrys->Write();
114  dPhivsThetaID->Write();
115  meanEnvsCrysID->Write();
116 
119  if (m_payloadName == "None") {
120  B2INFO("eclee5x5Algorithm has not been asked to find constants; copying input histograms and quitting");
121  histfile->Close();
122  return c_NotEnoughData;
123  }
124 
127  for (int cellID = 1; cellID <= 8736; cellID++) {
128 
130  if (AverageInitCalib->GetBinContent(cellID) > 0.) {
131  if (EntriesvsCrys->GetBinContent(cellID) < m_minEntries) {
132  histfile->Close();
133  B2INFO("eclee5x5Algorithm: insufficient data for cellID = " << cellID << " " << EntriesvsCrys->GetBinContent(cellID) << " entries");
134  return c_NotEnoughData;
135  }
136  }
137  }
138 
139 
143 
144 
147  bool foundConst = false;
148  TH1F* gVsCrysID = new TH1F("gVsCrysID", "Ratio of new to old calibration vs crystal ID;crystal ID;vector g", 8736, 0, 8736);
149  TH1F* CalibVsCrysID = new TH1F("CalibVsCrysID", "Calibration constant vs crystal ID;crystal ID;counts per GeV", 8736, 0, 8736);
150  TH1F* ExpEnergyperCrys = new TH1F("ExpEnergyperCrys", "Expected energy per crystal;Crystal ID;Energy in 25 crystal sum (GeV)", 8736,
151  0, 8736);
152  TString title = "dPhi cut per crystal " + m_payloadName + ";Crystal ID;dPhi requirement (deg)";
153  TH1F* dPhiperCrys = new TH1F("dPhiperCrys", title, 8736, 0, 8736);
154 
157  if (m_payloadName == "ECLExpee5x5E") {
158  for (int cellID = 1; cellID <= 8736; cellID++) {
159  float mean = meanEnvsCrysID->GetBinContent(cellID);
160  float stdDev = meanEnvsCrysID->GetBinError(cellID);
161  float inputE = AverageExpECrys->GetBinContent(cellID);
162  if (mean > 0. and stdDev > 0.) {
163  ExpEnergyperCrys->SetBinContent(cellID, mean * abs(inputE));
164  ExpEnergyperCrys->SetBinError(cellID, stdDev * abs(inputE));
165  } else {
166  ExpEnergyperCrys->SetBinContent(cellID, inputE);
167  ExpEnergyperCrys->SetBinError(cellID, 0.05 * inputE);
168  }
169  }
170 
171  //..Generate the payload, if requested
172  foundConst = true;
173  if (m_storeConst) {
174  std::vector<float> tempCalib;
175  std::vector<float> tempCalibStdDev;
176 
177  for (int cellID = 1; cellID <= 8736; cellID++) {
178  tempCalib.push_back(ExpEnergyperCrys->GetBinContent(cellID));
179  tempCalibStdDev.push_back(ExpEnergyperCrys->GetBinError(cellID));
180  }
181  ECLCrystalCalib* ExpectedE = new ECLCrystalCalib();
182  ExpectedE->setCalibVector(tempCalib, tempCalibStdDev);
183  saveCalibration(ExpectedE, "ECLExpee5x5E");
184  B2INFO("eclCosmicEAlgorithm: successfully stored expected energies ECLExpee5x5E");
185  }
186 
189  } else if (m_payloadName == "ECLCrystalEnergy5x5") {
190 
191  //..Create the Q matrix and the R vector
192  TMatrixDSym matrixQ(8736);
193  TVectorD vectorR(8736);
194  for (int ix = 1; ix <= 8736; ix++) {
195  vectorR[ix - 1] = RvsCrysID->GetBinContent(ix);
196  for (int iy = 1; iy <= 8736; iy++) {
197  matrixQ[ix - 1][iy - 1] = Qmatrix->GetBinContent(ix, iy);
198  }
199  }
200 
201  //..Crystals that are not being calibrated have no entries in NR. Adjust R and Q to get g=-1
202  int nNotCalibrated = 0;
203  for (int cellID = 1; cellID <= 8736; cellID++) {
204  if (NRvsCrysID->GetBinContent(cellID) < 0.5) {
205  for (int othercell = 1; othercell <= 8736; othercell++) {
206  matrixQ[cellID - 1][othercell - 1] = 0.;
207  matrixQ[othercell - 1][cellID - 1] = 0.;
208  }
209  matrixQ[cellID - 1][cellID - 1] = 1.;
210  vectorR[cellID - 1] = -1.;
211  nNotCalibrated++;
212  }
213  }
214  B2INFO("eclCosmicEAlgorithm: " << nNotCalibrated << " crystals will not be calibrated. ");
215 
216  //..Invert to solve Q g = R
217  TDecompLU lu(matrixQ);
218  bool solved;
219  TVectorD vectorg = lu.Solve(vectorR, solved);
220 
221  //..Fill histograms and check that there are no unexpected negative output values
222  if (solved) {
223  foundConst = true;
224  for (int cellID = 1; cellID <= 8736; cellID++) {
225  gVsCrysID->SetBinContent(cellID, vectorg[cellID - 1]);
226  gVsCrysID->SetBinError(cellID, 0.);
227  float newCalib = vectorg[cellID - 1] * abs(AverageInitCalib->GetBinContent(cellID));
228  CalibVsCrysID->SetBinContent(cellID, newCalib);
229  CalibVsCrysID->SetBinError(cellID, 0.);
230 
231  if (vectorg[cellID - 1] < 0. and NRvsCrysID->GetBinContent(cellID) > 0.) {foundConst = false;}
232  }
233  }
234 
235  //..Generate the payload if requested, and if the matrix inversion worked
236  if (m_storeConst and foundConst) {
237  std::vector<float> tempCalib;
238  std::vector<float> tempCalibStdDev;
239 
240  for (int cellID = 1; cellID <= 8736; cellID++) {
241  tempCalib.push_back(CalibVsCrysID->GetBinContent(cellID));
242  tempCalibStdDev.push_back(CalibVsCrysID->GetBinError(cellID));
243  }
244  ECLCrystalCalib* e5x5ECalib = new ECLCrystalCalib();
245  e5x5ECalib->setCalibVector(tempCalib, tempCalibStdDev);
246  saveCalibration(e5x5ECalib, "ECLCrystalEnergyee5x5");
247  B2INFO("eclCosmicEAlgorithm: successfully stored calibration ECLCrystalEnergyee5x5");
248  }
249 
252  } else if (m_payloadName == "ECLeedPhiData" or m_payloadName == "ECLeedPhiMC") {
253 
254  //..Find mean and sigma of Gaussian fit to ThetaID projections with sufficient statistics
255  float dPhiCenter[69] = {};
256  float dPhiHalfWidth[69] = {};
257  int firstID = 0;
258  int lastID = 0;
259  const int nbins = dPhivsThetaID->GetNbinsX();
260  for (int ib = 1; ib <= nbins; ib++) {
261  TH1D* proj = dPhivsThetaID->ProjectionY("proj", ib, ib);
262  if (proj->Integral() < 100) {continue;}
263  int thetaID = (int)dPhivsThetaID->GetXaxis()->GetBinCenter(ib);
264  if (firstID == 0) { firstID = thetaID; }
265  lastID = thetaID;
266 
267  //..Find the fit limits
268  double fracHi = m_fracHiSym;
269  double nsigHi = m_nsigHiSym;
270  if (thetaID <= m_lastLoThetaID) {
271  fracHi = m_fracHiASym;
272  nsigHi = m_nsigHiASym;
273  }
274 
275  //..Fit range includes all bins with entries>m_fracLo*peak on the low side of the
276  // peak, and entries>m_fracHi*peak on the high side.
277  // i.e. low edge of bin on low side, high edge on high side = low edge of bin+1
278  double peak = proj->GetMaximum();
279  int nPhiBins = proj->GetNbinsX();
280  int binLo = 1;
281  do {
282  binLo++;
283  } while (proj->GetBinContent(binLo) < m_fracLo * peak and binLo < nPhiBins);
284  double xfitLo = proj->GetBinLowEdge(binLo);
285 
286  //..Start search for the upper edge of the fit range at 175 deg to avoid gamma gamma
287  // and e gamma peaks
288  int binHi = proj->GetXaxis()->FindBin(175.01);
289  do {
290  binHi--;
291  } while (proj->GetBinContent(binHi)<fracHi* peak and binHi>1);
292  double xfitHi = proj->GetBinLowEdge(binHi + 1);
293 
294  //..Check that the fit region is sensible
295  int peakBin = proj->GetMaximumBin();
296  if (binLo >= peakBin or binHi <= peakBin) {
297  B2ERROR("Flawed dPhi fit range for thetaID = " << thetaID << " peakBin = " << peakBin << " binLo = " << binLo << "binHi = " <<
298  binHi);
299  }
300 
301  //..Now fit a Gaussian to the selected region
302  proj->Fit("gaus", "", "", xfitLo, xfitHi);
303 
304  //..Find mean, sigma, and selection range. Record center and half-width of range.
305  TF1* fitGaus = proj->GetFunction("gaus");
306  float mean = fitGaus->GetParameter(1);
307  float sigma = fitGaus->GetParameter(2);
308  float dPhiLo = mean - m_nsigLo * sigma;
309  float dPhiHi = mean + nsigHi * sigma;
310  dPhiCenter[thetaID] = 0.5 * (dPhiLo + dPhiHi);
311  dPhiHalfWidth[thetaID] = 0.5 * (dPhiHi - dPhiLo);
312  }
313 
314  //..Pad the thetaID's without fits with the first or last thetaID values
315  for (int thetaID = 0; thetaID < firstID; thetaID++) {
316  dPhiCenter[thetaID] = dPhiCenter[firstID];
317  dPhiHalfWidth[thetaID] = dPhiHalfWidth[firstID];
318  }
319  for (int thetaID = lastID + 1; thetaID < 69; thetaID++) {
320  dPhiCenter[thetaID] = dPhiCenter[lastID];
321  dPhiHalfWidth[thetaID] = dPhiHalfWidth[lastID];
322  }
323 
324  //..Now copy these to each crystal to generate the payload and fill the output histogram.
325  // We will use ECLNeighours to get the number of crystals in each theta ring
327  std::vector<float> tempCalib;
328  std::vector<float> tempCalibWidth;
329  tempCalib.resize(8736);
330  tempCalibWidth.resize(8736);
331  int crysID = 0;
332  for (int thetaID = 0; thetaID < 69; thetaID++) {
333  for (int ic = 0; ic < m_eclNeighbours5x5->getCrystalsPerRing(thetaID); ic++) {
334  tempCalib.at(crysID) = dPhiCenter[thetaID];
335  tempCalibWidth.at(crysID) = dPhiHalfWidth[thetaID];
336  dPhiperCrys->SetBinContent(crysID + 1, dPhiCenter[thetaID]);
337  dPhiperCrys->SetBinError(crysID + 1, dPhiHalfWidth[thetaID]);
338  crysID++;
339  }
340  }
341 
342  //..Store the payload, if requested
343  foundConst = true;
344  if (m_storeConst) {
345  ECLCrystalCalib* eedPhi = new ECLCrystalCalib();
346  eedPhi->setCalibVector(tempCalib, tempCalibWidth);
348  B2INFO("eclCosmicEAlgorithm: successfully stored calibration " << m_payloadName);
349  }
350 
353  } else {
354  B2ERROR("eclee5x5Algorithm: invalid payload name: m_payloadName = " << m_payloadName);
355  }
356 
359  histfile->cd();
360  if (m_payloadName == "ECLExpee5x5E") {
361  ExpEnergyperCrys->Write();
362  } else if (m_payloadName == "ECLCrystalEnergy5x5") {
363  gVsCrysID->Write();
364  CalibVsCrysID->Write();
365  } else if (m_payloadName == "ECLeedPhiData" or m_payloadName == "ECLeedPhiMC") {
366  dPhiperCrys->Write();
367  }
368  histfile->Close();
369 
370  dummy = (TH1F*)gROOT->FindObject("gVsCrysID"); delete dummy;
371  dummy = (TH1F*)gROOT->FindObject("CalibVsCrysID"); delete dummy;
372  dummy = (TH1F*)gROOT->FindObject("ExpEnergyperCrys"); delete dummy;
373  dummy = (TH1F*)gROOT->FindObject("dPhiperCrys"); delete dummy;
374 
377  if (!m_storeConst) {
378  if (foundConst) {
379  B2INFO("eclee5x5Algorithm successfully found constants but was not asked to store them");
380  } else {
381  B2INFO("eclee5x5Algorithm was not asked to store constants, and did not succeed in finding them");
382  }
383  return c_Failure;
384  } else if (!foundConst) {
385  if (m_payloadName == "ECLExpee5x5E") {
386  B2INFO("eclee5x5Algorithm: failed to store expected values");
387  } else if (m_payloadName == "ECLCrystalEnergy5x5") {
388  B2INFO("eclee5x5Algorithm: failed to store calibration constants");
389  } else if (m_payloadName == "ECLeedPhiData" or m_payloadName == "ECLeedPhiMC") {
390  B2INFO("eclee5x5Algorithm: failed to find dPhi* selection criteria");
391  }
392  return c_Failure;
393  }
394  return c_OK;
395 }
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.
@ c_NotEnoughData
Needs more data =2 in Python.
@ c_Failure
Failed =3 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.
Class to get the neighbours for a given cell id.
Definition: ECLNeighbours.h:23
short int getCrystalsPerRing(const short int thetaid) const
return number of crystals in a given theta ring
Definition: ECLNeighbours.h:37
int m_minEntries
all crystals to be calibrated must have this many entries
double m_fracHiASym
or fracHiASym*peak, at low values of thetaID
double m_fracLo
start dPhi fit where data is > fraclo*peak
double m_nsigHiASym
or mean+nsigHiASym*sigma at low thetaID
double m_nsigLo
dPhi region is mean - nsigLo*sigma
std::string m_payloadName
Name of the payload to be stored.
ECL::ECLNeighbours * m_eclNeighbours5x5
Neighbours, used to get nCrys per ring.
std::string m_outputName
..Parameters to control job to find energy calibration using Bhabhas
bool m_storeConst
write payload to localdb if true
int m_lastLoThetaID
use asymmetric dPhi range for thetaID<= this value
virtual EResult calibrate() override
..Run algorithm on events
double m_fracHiSym
end dPhi fit where data is > fracHiSym*peak
double m_nsigHiSym
to mean + nsigHiSym*sigma
Abstract base class for different kinds of events.