Belle II Software  release-08-00-10
eclLeakageAlgorithm.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/eclLeakageAlgorithm.h>
9 #include <ecl/calibration/tools.h>
10 #include <ecl/dbobjects/ECLLeakageCorrections.h>
11 #include <framework/datastore/StoreObjPtr.h>
12 #include <framework/dataobjects/EventMetaData.h>
13 #include <framework/datastore/DataStore.h>
14 
15 
16 #include "TH1D.h"
17 #include "TF1.h"
18 #include "TTree.h"
19 #include "TFile.h"
20 #include "TDirectory.h"
21 #include "TProfile.h"
22 #include <iostream>
23 
24 using namespace std;
25 using namespace Belle2;
26 using namespace ECL;
27 using namespace Calibration;
28 
29 /**************************************************************************
30  * eclLeakageAlgorithm analyzes single photon MC to find the eclLeakage payload
31  * After initial setup, the position dependent corrections are found using
32  * photons with good reconstruction.
33  * The correction that depends on the number of crystals has been replaced
34  * by the ECLnOptimal payload.
35  *
36  * Major steps in the process:
37  * 1. Set up (read in parameters and tree, define histogram binning)
38  * 2. Fill and fit histograms of normalized reconstructed energy, one per
39  * thetaID/energy. Novosibirsk fit parameter eta is floated in this fit, then
40  * fixed for fits of data with finer location binning. Peak (i.e. overall correction
41  * for this thetaID/energy) is used to normalize position-dependent correction.
42  * 3. Fill and fit histograms of normalized energy for each location (nominally 29
43  * locations in theta and phi per thetaID) to get the position-dependent correction.
44  * Nominally 69 x 8 x 29 x 3 = 48,024 histograms. 3 = theta, phi next to mech, phi not
45  * next to mech.
46  * 4. Pack payload quantities into vectors and histograms
47  * 5. Fill and fit summary and resolution histograms
48  * 6. Finish; close histogram file and and write payload
49 
50  **************************************************************************/
51 
52 
54 //..Function to get starting parameters for the Novo fit
55 std::vector<double> eclLeakageFitParameters(TH1F* h, const double& target)
56 {
57 
58  //..Find the smallest range of bins that contain at least "target" entries
59  double maxIntegral = h->Integral();
60  const int nBins = h->GetNbinsX();
61  int minLo = 1;
62  int minHi = nBins;
63 
64  //..Store a vector of integrals over histogram ranges so that I only need to
65  // call h->Integral nBins times, instead of nBins^2.
66  std::vector<double> intVector; // intVector[n] = sum of histogram bins [1,n+1] inclusive
67  intVector.push_back(h->GetBinContent(1));
68  for (int iLo = 2; iLo <= nBins; iLo++) {
69  double nextIntegral = intVector[iLo - 2] + h->GetBinContent(iLo);
70  intVector.push_back(nextIntegral);
71  }
72 
73  //..Now just look through all possible bin ranges to find the best one
74  for (int iLo = 2; iLo <= nBins; iLo++) {
75  for (int iHi = iLo; iHi <= nBins; iHi++) {
76 
77  //..sum[iLo, iHi] = sum[1, iHi] - sum[1,iLo-1] = intVector[iHi-1] - intVector[iLo-2]
78  double integral = intVector[iHi - 1] - intVector[iLo - 2];
79 
80  //..same number of bins, pick the pair with more entries
81  if ((integral > target and (iHi - iLo) < (minHi - minLo)) or
82  (integral > target and (iHi - iLo) == (minHi - minLo) and integral > maxIntegral)
83  ) {
84  minLo = iLo;
85  minHi = iHi;
86  maxIntegral = integral;
87  }
88  }
89  }
90 
91  //..Fit parameters are derived (in part) from this range of bins
92  const double lowE = h->GetBinLowEdge(minLo);
93  const double highE = h->GetBinLowEdge(minHi + 1);
94  const double peak = 0.5 * (lowE + highE);
95  const double sigma = 0.4 * (highE - lowE);
96  const double eta = 0.4; // typical value
97  const int nfitBins = (1 + minHi - minLo);
98 
99  std::vector<double> parameters;
100  parameters.push_back(peak);
101  parameters.push_back(sigma);
102  parameters.push_back(eta);
103  parameters.push_back(lowE);
104  parameters.push_back(highE);
105  parameters.push_back(nfitBins);
106 
107  return parameters;
108 }
109 
111 //..Novosibirsk; H. Ikeda et al. / NIM A 441 (2000) 401-426
112 // cppcheck-suppress constParameter ; TF1 fit functions cannot have const parameters
113 double eclLeakageNovo(Double_t* x, Double_t* par)
114 {
115 
116  Double_t peak = par[1];
117  Double_t width = par[2];
118  Double_t sln4 = sqrt(log(4));
119  Double_t y = par[3] * sln4;
120  Double_t tail = -log(y + sqrt(1 + y * y)) / sln4;
121  Double_t qc = 0.;
122 
123  if (TMath::Abs(tail) < 1.e-7)
124  qc = 0.5 * TMath::Power(((x[0] - peak) / width), 2);
125  else {
126  double qa = tail * sqrt(log(4.));
127  double qb = sinh(qa) / qa;
128  double qx = (x[0] - peak) / width * qb;
129  double qy = 1. + tail * qx;
130 
131  if (qy > 1.E-7)
132  qc = 0.5 * (TMath::Power((log(qy) / tail), 2) + tail * tail);
133  else
134  qc = 15.0;
135  }
136 
137  return par[0] * exp(-qc);
138 }
139 
140 
142 //..Function to check for bad fits
143 // 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb, 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit
144 int eclLeakageFitQuality(const double& fitLo, const double& fitHi, const double& fitPeak, const double& fitSigma,
145  const double& fitdEta, const double& fitProb)
146 {
147  const double tolerance = 0.02;
148  const double redoFitProb = 1e-5;
149  const double minFitProb = 1e-9;
151  int fitStatus = 0;
152  if (fitProb < redoFitProb) {fitStatus = 1;}
153  if (fitProb < minFitProb) {fitStatus = 3;}
154  double tol = tolerance * (fitHi - fitLo);
155  if (fitPeak < (fitLo + tol) or fitPeak > (fitHi - tol)) {fitStatus = 4;}
156  if (fitSigma<tol or fitSigma>(fitHi - fitLo - tol)) {fitStatus = 5;}
157  double tolEta = tolerance;
158  if (fitdEta < tolEta) {fitStatus = 6;} // just use the tolerance for eta (note that it is eta - limit that is passed)
159  return fitStatus;
160 }
161 
163 eclLeakageAlgorithm::eclLeakageAlgorithm(): CalibrationAlgorithm("eclLeakageCollector")
164 {
166  "Generate payload ECLLeakageCorrection from single photon MC recorded by eclLeakageCollectorModule"
167  );
168 }
169 
170 
171 
174 {
175 
176  //====================================================================================
177  //====================================================================================
178  //..Step 1. Set up prior to first loop through data
179  B2INFO("starting eclLeakageAlgorithm");
180 
181  //-----------------------------------------------------------------------------------
182  //..Read in histograms created by the collector and fix normalization
183  auto inputParameters = getObjectPtr<TH1F>("inputParameters");
184  int lastBin = inputParameters->GetNbinsX();
185  double scale = inputParameters->GetBinContent(lastBin); // number of times inputParameters was filled
186  for (int ib = 1; ib < lastBin; ib++) {
187  double param = inputParameters->GetBinContent(ib);
188  inputParameters->SetBinContent(ib, param / scale);
189  inputParameters->SetBinError(ib, 0.);
190  }
191 
192  //..And write to disk.
193  TFile* histfile = new TFile("eclLeakageAlgorithm.root", "recreate");
194  inputParameters->Write();
195 
196  //..Parameters
197  const int nThetaID = 69;
198  const int nLeakReg = 3;
199  const int firstBarrelID = 13;
200  const int lastBarrelID = 58;
201  const int firstUsefulThID = 3;
202  const int lastUsefulThID = 66;
204  const int nPositions = (int)(inputParameters->GetBinContent(1) + 0.001);
205  const int nEnergies = (int)(inputParameters->GetBinContent(2) + 0.001);
206  const int nbinX = nEnergies * nThetaID;
208  const double etaMin = -5.;
209  const double etaMax = 5.;
212  //..Energies
213  auto generatedE = create2Dvector<float>(nLeakReg, nEnergies); // 3 regions forward barrel backward
214 
215  int bin = 2; // bin 1 = nPositions, bin 2 = nEnergies, bin 3 = first energy
216  for (int ireg = 0; ireg < nLeakReg; ireg++) {
217  B2INFO("Generated energies for ireg = " << ireg);
218  for (int ie = 0; ie < nEnergies; ie++) {
219  bin++;
220  generatedE[ireg][ie] = inputParameters->GetBinContent(bin);
221  B2INFO(" " << ie << " " << generatedE[ireg][ie] << " GeV");
222  }
223  }
224  B2INFO("Low energy threshold = " << m_lowEnergyThreshold);
225 
226  //..Energy per thetaID (in MeV, for use in titles etc)
227  auto iEnergiesMeV = create2Dvector<int>(nEnergies, nThetaID);
228  for (int thID = 0; thID < nThetaID; thID++) {
229  int ireg = 0;
230  if (thID >= firstBarrelID and thID <= lastBarrelID) {
231  ireg = 1;
232  } else if (thID > lastBarrelID) {
233  ireg = 2;
234  }
235  for (int ie = 0; ie < nEnergies; ie++) {
236  iEnergiesMeV[ie][thID] = (int)(1000.*generatedE[ireg][ie] + 0.5);
237  }
238  }
239 
240  //-----------------------------------------------------------------------------------
241  //..Bins for eFrac histograms (eFrac = uncorrected reconstructed E / E true)
242  const double eFracLo = 0.4; // low edge of eFrac histograms
243  const double eFracHi = 1.5; // high edge of eFrac histograms
244  auto nEfracBins = create2Dvector<int>(nEnergies, nThetaID);
245  for (int thID = 0; thID < nThetaID; thID++) {
246  B2DEBUG(25, "eFrac nBins for thetaID " << thID);
247  for (int ie = 0; ie < nEnergies; ie++) {
248 
249  //..ballpark resolution
250  double res_squared = 0.0001 + 0.064 / iEnergiesMeV[ie][thID];
251 
252  //..Convert this to an even integer to get number of bins
253  double binNumOver2 = 3. / sqrt(res_squared);
254  int tempNBin = (int)(binNumOver2 + 0.5);
255  nEfracBins[ie][thID] = 2 * tempNBin;
256  B2DEBUG(25, " ie = " << ie << " E = " << iEnergiesMeV[ie][thID] << " nBins = " << nEfracBins[ie][thID]);
257  }
258  }
259 
260  //-----------------------------------------------------------------------------------
261  //..Set up Novosibirsk fit function
262  TString statusString[7] = {"good", "refit", "lowStat", "lowProb", "peakAtLimit", "sigmaAtLimit", "etaAtLimit"};
263  const double fracEnt[2] = {0.683, 0.5}; // fit range includes 68% or 50% of entries
264  const double minEntries = 100.; // don't use fits with fewer entries
265  const double minMaxBin = 50.; // rebin if max bin is below this value
266  const double highMaxBin = 300.; // can float eta if max bin is above this value
267 
268  TF1* func = new TF1("eclLeakageNovo", eclLeakageNovo, eFracLo, eFracHi, 4);
269  func->SetParNames("normalization", "peak", "effSigma", "eta");
270  func->SetLineColor(kRed);
271 
272  //-----------------------------------------------------------------------------------
273  //..Specify the TTree
274  auto tree = getObjectPtr<TTree>("tree");
275  tree->SetBranchAddress("cellID", &t_cellID);
276  tree->SetBranchAddress("thetaID", &t_thetaID);
277  tree->SetBranchAddress("region", &t_region);
278  tree->SetBranchAddress("thetaBin", &t_thetaBin);
279  tree->SetBranchAddress("phiBin", &t_phiBin);
280  tree->SetBranchAddress("phiMech", &t_phiMech);
281  tree->SetBranchAddress("energyBin", &t_energyBin);
282  tree->SetBranchAddress("nCrys", &t_nCrys);
283  tree->SetBranchAddress("energyFrac", &t_energyFrac);
284  tree->SetBranchAddress("origEnergyFrac", &t_origEnergyFrac);
285  tree->SetBranchAddress("locationError", &t_locationError);
286 
287  const int treeEntries = tree->GetEntries();
288  B2INFO("eclLeakageAlgorithm entries = " << treeEntries);
289 
290  //-----------------------------------------------------------------------------------
291  //..Get current payload to help validate the new payload
292 
293  //..Set event, run, exp number
294  const auto exprun = getRunList();
295  const int iEvt = 1;
296  const int iRun = exprun[0].second;
297  const int iExp = exprun[0].first;
300  evtPtr.registerInDataStore();
302  evtPtr.construct(iEvt, iRun, iExp);
303  DBStore& dbstore = DBStore::Instance();
304  dbstore.update();
305 
306  //..Existing payload
307  DBObjPtr<ECLLeakageCorrections> existingCorrections("ECLLeakageCorrections");
308  TH2F existingThetaCorrection = existingCorrections->getThetaCorrections();
309  existingThetaCorrection.SetName("existingThetaCorrection");
310  TH2F existingPhiCorrection = existingCorrections->getPhiCorrections();
311  existingPhiCorrection.SetName("existingPhiCorrection");
312 
313  //..Write out the correction histograms
314  histfile->cd();
315  existingThetaCorrection.Write();
316  existingPhiCorrection.Write();
317 
318 
319  //====================================================================================
320  //====================================================================================
321  //..Step 2. First loop, fill histograms of e/eTrue for each energy and thetaID
322 
323  //-----------------------------------------------------------------------------------
324  //..One histogram of e/eTrue per energy per thetaID.
325  // Used to fix eta in subsequent fits, and to get overall correction for that thetaID
326  // (which should be very close to 1).
327  TString name;
328  TString title;
329  auto hELabUncorr = create2Dvector<TH1F*>(nEnergies, nThetaID);
330  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
331  TString sthID = std::to_string(thID);
332  for (int ie = 0; ie < nEnergies; ie++) {
333  name = "hELabUncorr_" + std::to_string(ie) + "_" + sthID;
334  title = "eFrac " + to_string(iEnergiesMeV[ie][thID]) + " MeV thetaID " + sthID + ";E/Etrue";
335 
336  //..High statistics for these plots; use more bins
337  hELabUncorr[ie][thID] = new TH1F(name, title, 2 * nEfracBins[ie][thID], eFracLo, eFracHi);
338  }
339  }
340 
341  //-----------------------------------------------------------------------------------
342  //..Loop over events and store eFrac
343  for (int i = 0; i < treeEntries; i++) {
344  tree->GetEntry(i);
345 
346  //..Only events with good reconstruction
347  bool goodReco = t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0;
348  if (not goodReco) {continue;}
349 
350  //..Fill histogram for full thetaID
351  hELabUncorr[t_energyBin][t_thetaID]->Fill(t_energyFrac);
352  }
353 
354  //-----------------------------------------------------------------------------------
355  //..Fit each thetaID/energy histogram to get peak (overall correction) and eta (fixed
356  // in subsequent fits to individual locations).
357  auto peakUncorr = create2Dvector<float>(nEnergies, nThetaID); // store peak from each fit
358  auto etaUncorr = create2Dvector<float>(nEnergies, nThetaID); // store eta from each fit
359 
360  std::vector<TString> failedELabUncorr; // names of hists with failed fits
361  std::vector<int> statusELabUncorr; // status of failed fits
362  int payloadStatus = 0; // Overall status of payload determination
363 
364  //..Record fit status
365  TH1F* statusOfhELabUncorr = new TH1F("statusOfhELabUncorr",
366  "status of hELabUncorr fits for each thetaID/energy. 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb, 4 = peakAtLimit, 5 = sigmaAtLimit",
367  6, 0,
368  6);
369 
370 
371  //..Fit each histogram
372  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
373  TString sthID = std::to_string(thID);
374  for (int ie = 0; ie < nEnergies; ie++) {
375  TH1F* hEnergy = (TH1F*)hELabUncorr[ie][thID];
376  double peak = -1.;
377  double eta = 0.;
378  int fitStatus = 2;
379  double entries = hEnergy->Integral();
380  int nIter = 0; // keep track of attempts to fit this histogram
381  double genE = iEnergiesMeV[ie][thID] / 1000.;
382  bool fitHist = entries > minEntries;
383 
384  //..Possibly iterate fit starting from this point
385  while (fitHist) {
386 
387  //..Fit parameters
388  double norm = hEnergy->GetMaximum();
389  double target = fracEnt[nIter] * entries; // fit range contains 68% or 50%
390  std::vector<double> startingParameters;// peak, sigma, eta, fitLow, fitHigh, nbins
391  startingParameters = eclLeakageFitParameters(hEnergy, target);
392  func->SetParameters(norm, startingParameters[0], startingParameters[1], startingParameters[2]);
393  func->SetParLimits(1, startingParameters[3], startingParameters[4]);
394  func->SetParLimits(2, 0., startingParameters[4] - startingParameters[3]);
395  func->SetParLimits(3, etaMin, etaMax);
396 
397  //..Fit
398  name = hEnergy->GetName();
399  B2DEBUG(40, "Fitting " << name.Data());
400  hEnergy->Fit(func, "LIEQ", "", startingParameters[3], startingParameters[4]);
401  peak = func->GetParameter(1);
402  double effSigma = func->GetParameter(2);
403  eta = func->GetParameter(3);
404  double prob = func->GetProb();
405 
406  //..Check fit quality 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb,
407  // 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit
408  double dEta = min((etaMax - eta), (eta - etaMin));
409  fitStatus = eclLeakageFitQuality(startingParameters[3], startingParameters[4], peak, effSigma, dEta, prob);
410 
411  //..If the fit probability is low, refit using a smaller range (fracEnt)
412  if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
413  nIter++;
414  } else {
415  fitHist = false;
416  }
417  }
418 
419  //-----------------------------------------------------------------------------------
420  //..Record failures and fit results.
421  // Mark payload as failed if energy is above low-energy threshold.
422  statusOfhELabUncorr->Fill(fitStatus + 0.000001);
423  if (fitStatus == 2 or fitStatus >= 4) {
424  statusELabUncorr.push_back(fitStatus);
425  failedELabUncorr.push_back(hEnergy->GetName());
426  if (genE > m_lowEnergyThreshold) {
427  payloadStatus = 1; // algorithm failed
428  } else {
429  peak = -1.; // fix up later
430  }
431  }
432  peakUncorr[ie][thID] = peak;
433  etaUncorr[ie][thID] = eta;
434 
435  //..Write to disk
436  if (entries > minEntries) {
437  histfile->cd();
438  hELabUncorr[ie][thID]->Write();
439  }
440  }
441  }
442 
443  //..Write out summary of fit status
444  histfile->cd();
445  statusOfhELabUncorr->Write();
446 
447  //-----------------------------------------------------------------------------------
448  //..Quit now if one of the high-energy fits failed, since we will not get a payload.
449  int nbadFit = statusELabUncorr.size();
450  if (nbadFit > 0) {B2ERROR("hELabUncorr fit failures (one histogram per energy/thetaID):");}
451  for (int ibad = 0; ibad < nbadFit; ibad++) {
452  int badStat = statusELabUncorr[ibad];
453  B2ERROR(" histogram " << failedELabUncorr[ibad].Data() << " status " << badStat << " " << statusString[badStat].Data());
454  }
455  if (payloadStatus != 0) {
456  B2ERROR("ecLeakageAlgorithm: fit to hELabUncorr failed. ");
457  histfile->Close();
458  return c_Failure;
459  }
460 
461  //-----------------------------------------------------------------------------------
462  //..Fix up any failed (low-energy) fits by using a neighbouring thetaID
463  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
464  for (int ie = 0; ie < nEnergies; ie++) {
465  if (peakUncorr[ie][thID] < 0.) {
466  if (thID > 40) {
467  for (int nextID = thID - 1; thID >= firstUsefulThID; thID--) {
468  if (peakUncorr[ie][nextID] > 0.) {
469  peakUncorr[ie][thID] = peakUncorr[ie][nextID];
470  break;
471  }
472  }
473  } else {
474  for (int nextID = thID + 1; thID <= lastUsefulThID; thID++) {
475  if (peakUncorr[ie][nextID] > 0.) {
476  peakUncorr[ie][thID] = peakUncorr[ie][nextID];
477  break;
478  }
479  }
480  }
481  }
482 
483  //..If we were unable to get a successful fit from a neighbour, give up
484  if (peakUncorr[ie][thID] < 0.) {
485  B2ERROR("ecLeakageAlgorithm: unable to correct hELabUncorr for thetaID " << thID << " ie " << ie);
486  histfile->Close();
487  return c_Failure;
488  }
489  }
490  }
491 
492  //====================================================================================
493  //====================================================================================
494  //..Step 3. Second loop, fill histograms of e/eTrue as a function of position.
495  // 29 locations in theta, 29 in phi.
496  // Crystals next to mechanical structure in phi are treated separately from
497  // crystals without mechanical structure.
498 
499  //-----------------------------------------------------------------------------------
500  //..Histograms to store the energy
501  const int nDir = 3;
502  const TString dirName[nDir] = {"theta", "phiMech", "phiNoMech"};
503 
504  auto eFracPosition = create4Dvector<TH1F*>(nEnergies, nThetaID, nDir, nPositions); // the histograms
505 
506 
507  std::vector<TString> failedeFracPosition; // names of hists with failed fits
508  std::vector<int> statuseFracPosition; // status of failed fits
509 
510  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
511  TString sthID = std::to_string(thID);
512  for (int ie = 0; ie < nEnergies; ie++) {
513  TString sie = std::to_string(ie);
514  for (int idir = 0; idir < nDir; idir++) {
515  TString sidir = std::to_string(idir);
516  for (int ipos = 0; ipos < nPositions; ipos++) {
517  TString sipos = std::to_string(ipos);
518  name = "eFracPosition_" + sie + "_" + sthID + "_" + sidir + "_" + sipos;
519  title = "eFrac " + to_string(iEnergiesMeV[ie][thID]) + " MeV thetaID " + sthID + " " + dirName[idir] + " pos " + sipos +
520  "; E/Etrue";
521  eFracPosition[ie][thID][idir][ipos] = new TH1F(name, title, nEfracBins[ie][thID], eFracLo, eFracHi);
522  }
523  }
524  }
525  }
526 
527  //..And some summary histograms
528  TH1F* statusOfeFracPosition = new TH1F("statusOfeFracPosition",
529  "eFrac fit status: -2 noData, -1 lowE, 0 good, 1 redo fit, 2 lowStat, 3 lowProb, 4 peakAtLimit, 5 sigmaAtLimit", 8, -2,
530  6);
531  TH1F* probOfeFracPosition = new TH1F("probOfeFracPosition", "fit probability of eFrac fits for each position;probability", 100, 0,
532  1);
533  TH1F* maxOfeFracPosition = new TH1F("maxOfeFracPosition", "max entries of eFrac histograms;maximum bin content", 100, 0, 1000);
534 
535  //-----------------------------------------------------------------------------------
536  //..Loop over events and store eFrac
537  for (int i = 0; i < treeEntries; i++) {
538  tree->GetEntry(i);
539 
540  //..Only events with good reconstruction
541  bool goodReco = t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0;
542  if (not goodReco) {continue;}
543 
544  //..Theta location (idir = 0)
545  int idir = 0;
546  eFracPosition[t_energyBin][t_thetaID][idir][t_thetaBin]->Fill(t_energyFrac);
547 
548  //..phi location. Note that t_phiBin starts at mechanical edge, if there is one.
549  idir = t_phiMech + 1;
550  eFracPosition[t_energyBin][t_thetaID][idir][t_phiBin]->Fill(t_energyFrac);
551  }
552 
553  //-----------------------------------------------------------------------------------
554  //..Now fit many many histograms to get the position-dependent corrections
555  auto positionCorrection = create4Dvector<float>(nEnergies, nThetaID, nDir, nPositions);
556  auto positionCorrectionUnc = create4Dvector<float>(nEnergies, nThetaID, nDir, nPositions);
557  int nHistToFit = 0;
558 
559 
560  //..Temp histogram of position corrections
561  TH1F* thetaCorSummary = new TH1F("thetaCorSummary", "Theta dependent corrections;theta dependent correction", 100, 0.4, 1.4);
562  TH1F* phiCorSummary = new TH1F("phiCorSummary", "Phi dependent corrections;phi dependent correction", 100, 0.4, 1.4);
563 
564  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
565  for (int ie = 0; ie < nEnergies; ie++) {
566  double genE = iEnergiesMeV[ie][thID] / 1000.;
567  for (int idir = 0; idir < nDir; idir++) {
568  for (int ipos = 0; ipos < nPositions; ipos++) {
569  TH1F* hEnergy = (TH1F*)eFracPosition[ie][thID][idir][ipos];
570  if (hEnergy->Integral() > minEntries) {nHistToFit++;}
571 
572  //..Default peak from fit to full thetaID/energy = peakUncorr;
573  // correction = peak / sqrt(peakUncorr) = sqrt(peakUncorr)
574  double correction = sqrt(peakUncorr[ie][thID]);
575  double corrUnc = 0.05; // arbitrary uncertainty
576  double prob = -1.;
577  int fitStatus = 2; // low stats
578  if (genE < m_lowEnergyThreshold) {fitStatus = -1;} // low energy, don't fit, use default peak value
579  if (hEnergy->GetEntries() < 0.5) {fitStatus = -2;} // unused, eg barrel pos=2
580  int nIter = 0; // keep track of attempts to fit this histogram
581  double entries = hEnergy->Integral();
582  bool fitHist = entries > minEntries and genE >= m_lowEnergyThreshold;
583 
584  //..Possibly iterate fit starting from this point
585  while (fitHist) {
586  fitHist = false;
587 
588  //..Fit parameters
589  double target = fracEnt[nIter] * entries; // fit range contains 68%
590  std::vector<double> startingParameters;// peak, sigma, eta, fitLow, fitHigh, nbins
591  bool findParm = true;
592  while (findParm) {
593  startingParameters = eclLeakageFitParameters(hEnergy, target);
594 
595  //..Rebin if stats are low and we have enough DOF
596  if (hEnergy->GetMaximum()<minMaxBin and startingParameters[5]>11.9) {
597  hEnergy->Rebin(2);
598  } else {
599  findParm = false;
600  }
601  }
602  double norm = hEnergy->GetMaximum(); // max bin after rebinning
603  double fitLow = startingParameters[3];
604  double fitHigh = startingParameters[4];
605 
606  //..Eta from the fit to the full crystal
607  double etaFix = etaUncorr[ie][thID];
608  func->SetParameters(norm, startingParameters[0], startingParameters[1], etaFix);
609  func->SetParLimits(1, fitLow, fitHigh);
610  func->SetParLimits(2, 0., fitHigh - fitLow);
611 
612  //..Fix eta, unless really good statistics
613  if (hEnergy->GetMaximum() < highMaxBin) {
614  func->SetParLimits(3, etaFix, etaFix);
615  } else {
616  func->SetParLimits(3, etaMin, etaMax);
617  }
618 
619  //..Fit
620  name = hEnergy->GetName();
621  B2DEBUG(40, "Fitting " << name.Data());
622  hEnergy->Fit(func, "LIEQ", "", fitLow, fitHigh);
623  double peak = func->GetParameter(1);
624  double effSigma = func->GetParameter(2);
625  double eta = func->GetParameter(3);
626  prob = func->GetProb();
627 
628  //..Check fit quality 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb,
629  // 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit.
630  double dEta = min((etaMax - eta), (eta - etaMin));
631  fitStatus = eclLeakageFitQuality(fitLow, fitHigh, peak, effSigma, dEta, prob);
632 
633  //..If the fit probability is low, refit using a smaller range (fracEnt)
634  if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
635  nIter++;
636  fitHist = true;
637 
638  // Store correction except for peak or sigma at limit.
639  } else if (fitStatus <= 3) {
640 
641  //..Divide by sqrt(peak for full crystal) to get correct average
642  // when multiplying the theta and phi corrections
643  correction = peak / sqrt(peakUncorr[ie][thID]);
644  corrUnc = func->GetParError(1) / sqrt(peakUncorr[ie][thID]);
645  }
646  }
647 
648  //..Store the correction for this position
649  positionCorrection[ie][thID][idir][ipos] = correction;
650  positionCorrectionUnc[ie][thID][idir][ipos] = corrUnc;
651  statusOfeFracPosition->Fill(fitStatus + 0.00001);
652  probOfeFracPosition->Fill(prob);
653  maxOfeFracPosition->Fill(hEnergy->GetMaximum());
654 
655  //..Summary histogram for successful fits. Note that fitStatus<0 also has prob<0.
656  if (prob > 0 and fitStatus <= 3) {
657  if (idir == 0) {
658  thetaCorSummary->Fill(correction);
659  } else {
660  phiCorSummary->Fill(correction);
661  }
662  }
663 
664  //..Record the failed fit details. We may still use the correction.
665  if (fitStatus >= 2) {
666  statuseFracPosition.push_back(fitStatus);
667  failedeFracPosition.push_back(hEnergy->GetName());
668  histfile->cd();
669  hEnergy->Write();
670  }
671  }
672  }
673  }
674  }
675 
676  //-----------------------------------------------------------------------------------
677  //..Summarize position fit results
678  B2INFO("eclLeakageAlgorithm: " << nHistToFit << " eFracPosition histograms to fit");
679  nbadFit = statuseFracPosition.size();
680  B2INFO("eFracPosition failed fits: " << nbadFit);
681  for (int ibad = 0; ibad < nbadFit; ibad++) {
682  int badStat = statuseFracPosition[ibad];
683  B2ERROR(" histogram " << failedeFracPosition[ibad].Data() << " status " << badStat << " " << statusString[badStat].Data());
684  }
685 
686  //..Write to disk
687  histfile->cd();
688  statusOfeFracPosition->Write();
689  probOfeFracPosition->Write();
690  maxOfeFracPosition->Write();
691  thetaCorSummary->Write();
692  phiCorSummary->Write();
693 
694 
695  //====================================================================================
696  //====================================================================================
697  //..Step 4. Store quantities needed for the payloads
698 
699  //-----------------------------------------------------------------------------------
700  //..First, we need to fix up the thetaID that have been so far missed.
701  // Just copy the values from the first or last thetaID for which corrections were found.
702 
703  //..ThetaID before the first useful one
704  for (int thID = firstUsefulThID - 1; thID >= 0; thID--) {
705  for (int ie = 0; ie < nEnergies; ie++) {
706 
707  for (int idir = 0; idir < 3; idir++) {
708  for (int ipos = 0; ipos < nPositions; ipos++) {
709  positionCorrection[ie][thID][idir][ipos] = positionCorrection[ie][thID + 1][idir][ipos];
710  positionCorrectionUnc[ie][thID][idir][ipos] = positionCorrectionUnc[ie][thID + 1][idir][ipos];
711  }
712  }
713  }
714  }
715 
716  //..ThetaID beyond last useful one
717  for (int thID = lastUsefulThID + 1; thID < nThetaID; thID++) {
718  for (int ie = 0; ie < nEnergies; ie++) {
719 
720  for (int idir = 0; idir < 3; idir++) {
721  for (int ipos = 0; ipos < nPositions; ipos++) {
722  positionCorrection[ie][thID][idir][ipos] = positionCorrection[ie][thID - 1][idir][ipos];
723  positionCorrectionUnc[ie][thID][idir][ipos] = positionCorrectionUnc[ie][thID - 1][idir][ipos];
724  }
725  }
726  }
727  }
728 
729  //-----------------------------------------------------------------------------------
730  //..std::vectors of generated energies
731  std::vector<float> forwardVector;
732  std::vector<float> barrelVector;
733  std::vector<float> backwardVector;
734  for (int ie = 0; ie < nEnergies; ie++) {
735  forwardVector.push_back(log(generatedE[0][ie]));
736  barrelVector.push_back(log(generatedE[1][ie]));
737  backwardVector.push_back(log(generatedE[2][ie]));
738  }
739 
740  //..Store in array format for validation studies
741  float leakLogE[3][12] = {};
742  for (int ie = 0; ie < nEnergies; ie++) {
743  leakLogE[0][ie] = forwardVector[ie];
744  leakLogE[1][ie] = barrelVector[ie];
745  leakLogE[2][ie] = backwardVector[ie];
746  }
747 
748  //-----------------------------------------------------------------------------------
749  //..Position dependent corrections
750 
751  //..2D histogram of theta-dependent corrections
752  TH2F thetaCorrection("thetaCorrection", "Theta correction vs bin;bin = thetaID + 69*energyBin;local theta bin", nbinX, 0, nbinX,
753  nPositions, 0, nPositions);
754  int ix = 0;
755  for (int ie = 0; ie < nEnergies; ie++) {
756  for (int thID = 0; thID < nThetaID; thID++) {
757  ix++;
758  for (int ipos = 0; ipos < nPositions; ipos++) {
759  int iy = ipos + 1;
760  float correction = positionCorrection[ie][thID][0][ipos];
761  float corrUnc = positionCorrectionUnc[ie][thID][0][ipos];
762  thetaCorrection.SetBinContent(ix, iy, correction);
763  thetaCorrection.SetBinError(ix, iy, corrUnc);
764  }
765  }
766  }
767 
768  //..2D histogram of phi-dependent corrections. Twice as many x bins;
769  // one set for crystals next to mechanical structure, one for otherwise.
770  TH2F phiCorrection("phiCorrection", "Phi correction vs bin;bin = thetaID + 69*energyBin;local phi bin", 2 * nbinX, 0, 2 * nbinX,
771  nPositions, 0, nPositions);
772  ix = 0;
773  for (int ie = 0; ie < nEnergies; ie++) {
774  for (int thID = 0; thID < nThetaID; thID++) {
775  ix++;
776  for (int ipos = 0; ipos < nPositions; ipos++) {
777  int iy = ipos + 1;
778 
779  //..First set of corrections
780  float correction = positionCorrection[ie][thID][1][ipos];
781  float corrUnc = positionCorrectionUnc[ie][thID][1][ipos];
782  phiCorrection.SetBinContent(ix, iy, correction);
783  phiCorrection.SetBinError(ix, iy, corrUnc);
784 
785  //..Second set
786  correction = positionCorrection[ie][thID][2][ipos];
787  corrUnc = positionCorrectionUnc[ie][thID][2][ipos];
788  phiCorrection.SetBinContent(ix + nbinX, iy, correction);
789  phiCorrection.SetBinError(ix + nbinX, iy, corrUnc);
790  }
791  }
792  }
793 
794  //-----------------------------------------------------------------------------------
795  //..nCrys dependent corrections; not used since nOptimal was implemented.
796  const int maxN = 21; // maximum number of crystals in a cluster
797  TH2F nCrystalCorrection("nCrystalCorrection", "nCrys correction vs bin;bin = thetaID + 69*energyBin;nCrys", nbinX, 0, nbinX,
798  maxN + 1, 0, maxN + 1);
799 
800  ix = 0;
801  for (int ie = 0; ie < nEnergies; ie++) {
802  for (int thID = 0; thID < nThetaID; thID++) {
803  ix++;
804  for (int in = 0; in <= maxN; in++) {
805  int iy = in + 1;
806  float correction = 1.;
807  float corrUnc = 0.;
808  nCrystalCorrection.SetBinContent(ix, iy, correction);
809  nCrystalCorrection.SetBinError(ix, iy, corrUnc);
810  }
811  }
812  }
813 
814 
815  //====================================================================================
816  //====================================================================================
817  //..Step 5. Diagnostic histograms
818 
819  //..Start by storing the payload histograms
820  histfile->cd();
821  thetaCorrection.Write();
822  phiCorrection.Write();
823  nCrystalCorrection.Write();
824 
825  //-----------------------------------------------------------------------------------
826  //..One histogram of new and original reconstructed energy after leakage correction
827  // per generated energy per region. Also uncorrected.
828  // The "Corrected no nCrys" and "Corrected measured" are identical, but keep both
829  // for easier comparison with earlier results.
830  const int nResType = 5;
831  const TString resName[nResType] = {"Uncorrected", "Original", "Corrected no nCrys", "Corrected measured", "Corrected true"};
832  const TString regName[nLeakReg] = {"forward", "barrel", "backward"};
833  auto energyResolution = create3Dvector<TH1F*>(nLeakReg, nEnergies, nResType);
834 
835  //..Base number of bins on a typical thetaID for each region
836  int thIDReg[nLeakReg];
837  thIDReg[0] = 9;
838  thIDReg[1] = 35;
839  thIDReg[2] = 61;
840 
841  for (int ireg = 0; ireg < nLeakReg; ireg++) {
842  for (int ie = 0; ie < nEnergies; ie++) {
843 
844  //..Titles and bins for this region and energy
845  TString stireg = std::to_string(ireg);
846  TString stie = std::to_string(ie);
847  TString stieName = std::to_string(iEnergiesMeV[ie][thIDReg[ireg]]);
848  int nbinReg = 20 * nEfracBins[ie][thIDReg[ireg]];
849 
850  for (int ires = 0; ires < nResType; ires++) {
851  name = "energyResolution_" + stireg + "_" + stie + "_" + std::to_string(ires);
852  title = resName[ires] + " energy, " + stieName + " MeV, " + regName[ireg] + ";Reconstructed E/Etrue";
853  energyResolution[ireg][ie][ires] = new TH1F(name, title, nbinReg, eFracLo, eFracHi);
854  }
855  }
856  }
857 
858  //-----------------------------------------------------------------------------------
859  //..Loop over events and store leakage-corrected energies
860  for (int i = 0; i < treeEntries; i++) {
861  tree->GetEntry(i);
862 
863  //..Only events with good reconstruction
864  bool goodReco = t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0;
865  if (not goodReco) {continue;}
866 
867  //-----------------------------------------------------------------------------------
868  //..Corrections using true energy
869 
870  //..Position-dependent leakage corrections using true energy
871  int idir = 0;
872  float thetaPosCor = positionCorrection[t_energyBin][t_thetaID][idir][t_thetaBin];
873 
874  idir = t_phiMech + 1;
875  float phiPosCor = positionCorrection[t_energyBin][t_thetaID][idir][t_phiBin];
876  float corrTrue = thetaPosCor * phiPosCor;
877 
878  //-----------------------------------------------------------------------------------
879  //..Find correction using measured energy. The correction is a function of corrected
880  // energy, so will need to iterate
881  float corrMeasured = 0.96; // typical correction as starting point
882  for (int iter = 0; iter < 2; iter++) {
883 
884 
885  //..Energy points that bracket this value
886  float energyRaw = t_energyFrac * generatedE[t_region][t_energyBin];
887  float logEnergy = log(energyRaw / corrMeasured);
888  int ie0 = 0; // lower energy point
889  if (logEnergy < leakLogE[t_region][0]) {
890  ie0 = 0;
891  } else if (logEnergy > leakLogE[t_region][nEnergies - 1]) {
892  ie0 = nEnergies - 2;
893  } else {
894  while (logEnergy > leakLogE[t_region][ie0 + 1]) {ie0++;}
895  }
896 
897  //..Corrections from lower and upper energy points
898  float cor0 = positionCorrection[ie0][t_thetaID][0][t_thetaBin] * positionCorrection[ie0][t_thetaID][t_phiMech + 1][t_phiBin];
899  float cor1 = positionCorrection[ie0 + 1][t_thetaID][0][t_thetaBin] * positionCorrection[ie0 + 1][t_thetaID][t_phiMech +
900  1][t_phiBin];
901 
902  //..Interpolate in logE
903  corrMeasured = cor0 + (cor1 - cor0) * (logEnergy - leakLogE[t_region][ie0]) / (leakLogE[t_region][ie0 + 1] -
904  leakLogE[t_region][ie0]);
905  }
906 
907  //..No longer have a separate case that excludes the nCrys correction
908  float corrNonCrys = corrMeasured;
909 
910  //-----------------------------------------------------------------------------------
911  //..Fill the histograms
912  energyResolution[t_region][t_energyBin][0]->Fill(t_energyFrac); // uncorrected
913  energyResolution[t_region][t_energyBin][1]->Fill(t_origEnergyFrac); // original payload
914  energyResolution[t_region][t_energyBin][2]->Fill(t_energyFrac / corrNonCrys); // no nCrys, measured energy
915  energyResolution[t_region][t_energyBin][3]->Fill(t_energyFrac / corrMeasured); // corrected, measured energy
916  energyResolution[t_region][t_energyBin][4]->Fill(t_energyFrac / corrTrue); // corrected, true energy
917  }
918 
919  //-----------------------------------------------------------------------------------
920  //..Fit each histogram to find peak, and extract resolution.
921 
922  //..Store the peak and resolution values for each histogram
923  auto peakEnergy = create3Dvector<float>(nLeakReg, nEnergies, nResType);
924  auto energyRes = create3Dvector<float>(nLeakReg, nEnergies, nResType);
925 
926  //..Loop over the histograms to be fit
927  for (int ireg = 0; ireg < nLeakReg; ireg++) {
928  for (int ie = 0; ie < nEnergies; ie++) {
929 
930  //..Base the resolution on the number of entries in the corrected plot for all 4 cases.
931  double entries = energyResolution[ireg][ie][nResType - 1]->Integral();
932  for (int ires = 0; ires < nResType; ires++) {
933  TH1F* hEnergy = (TH1F*)energyResolution[ireg][ie][ires];
934  double peak = -1.;
935  double res68 = 99.; // resolution based on 68.3% of the entries
936  int nIter = 0; // keep track of attempts to fit this histogram
937  bool fitHist = entries > minEntries;
938 
939  //..Possibly iterate fit starting from this point
940  while (fitHist) {
941 
942  //..Fit parameters
943  double norm = hEnergy->GetMaximum();
944  double target = fracEnt[nIter] * entries; // fit range contains 68.3% or 50%
945  std::vector<double> startingParameters;// peak, sigma, eta, fitLow, fitHigh, nbins
946  startingParameters = eclLeakageFitParameters(hEnergy, target);
947  func->SetParameters(norm, startingParameters[0], startingParameters[1], startingParameters[2]);
948  func->SetParLimits(1, startingParameters[3], startingParameters[4]);
949  func->SetParLimits(2, 0., startingParameters[4] - startingParameters[3]);
950  func->SetParLimits(3, etaMin, etaMax);
951 
952  //..First iteration the fitting range contains >68.3% of the histogram integral.
953  if (nIter == 0) {
954 
955  //..Find the bin numbers used in the fit
956  int minLo = hEnergy->GetXaxis()->FindBin(startingParameters[3] + 0.0001);
957  int minHi = hEnergy->GetXaxis()->FindBin(startingParameters[4] - 0.0001);
958 
959  //..Subtract a partial bin to get to exactly 68.3%
960  double dx = hEnergy->GetBinLowEdge(minLo + 1) - hEnergy->GetBinLowEdge(minLo);
961  double overage = hEnergy->Integral(minLo, minHi) - target;
962  double subLo = overage / hEnergy->GetBinContent(minLo);
963  double subHi = overage / hEnergy->GetBinContent(minHi);
964 
965  //..Resolution is half the range that contains 68.3% of the events
966  res68 = 0.5 * dx * (1 + minHi - minLo - max(subLo, subHi));
967  }
968 
969  //..Fit
970  name = hEnergy->GetName();
971  B2DEBUG(40, "Fitting " << name.Data());
972  hEnergy->Fit(func, "LIEQ", "", startingParameters[3], startingParameters[4]);
973  peak = func->GetParameter(1);
974  double effSigma = func->GetParameter(2);
975  double eta = func->GetParameter(3);
976  double prob = func->GetProb();
977 
978  //..Check fit quality 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb,
979  // 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit
980  double dEta = min((etaMax - eta), (eta - etaMin));
981  int fitStatus = eclLeakageFitQuality(startingParameters[3], startingParameters[4], peak, effSigma, dEta, prob);
982 
983  //..If the fit probability is low, refit using a smaller range (fracEnt)
984  if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
985  nIter++;
986  } else {
987  fitHist = false;
988  }
989 
990  } // end of while
991 
992  //..Record peak and resolution
993  peakEnergy[ireg][ie][ires] = peak;
994  energyRes[ireg][ie][ires] = res68;
995 
996  //..And write to disk
997  histfile->cd();
998  hEnergy->Write();
999  }
1000  }
1001  }
1002 
1003  //-----------------------------------------------------------------------------------
1004  //..Summarize resolution
1005  int nresBins = nEnergies * nLeakReg * (nResType + 1); // +1 to add an empty bin after each set
1006  TH1F* peakSummary = new TH1F("peakSummary", "Peak E/Etrue for each method, region, energy;Energy energy point", nresBins, 0,
1007  nEnergies);
1008  TH1F* resolutionSummary = new TH1F("resolutionSummary", "Resolution/peak for each method, region, energy;Energy energy point",
1009  nresBins, 0, nEnergies);
1010 
1011  B2INFO("Resolution divided by peak for each energy bin and region " << nResType << " ways");
1012  for (int ires = 0; ires < nResType; ires++) {B2INFO(" " << resName[ires]);}
1013  ix = 0;
1014  for (int ie = 0; ie < nEnergies; ie++) {
1015  B2INFO(" energy point " << ie);
1016  for (int ireg = 0; ireg < nLeakReg; ireg++) {
1017  B2INFO(" region " << ireg);
1018  for (int ires = 0; ires < nResType; ires++) {
1019 
1020  //..Store in summary histograms
1021  ix++;
1022  peakSummary->SetBinContent(ix, peakEnergy[ireg][ie][ires]);
1023  peakSummary->SetBinError(ix, 0.);
1024  resolutionSummary->SetBinContent(ix, energyRes[ireg][ie][ires] / peakEnergy[ireg][ie][ires]);
1025  resolutionSummary->SetBinError(ix, 0.);
1026 
1027  //..Print out as well
1028  B2INFO(" " << ires << " " << energyRes[ireg][ie][ires] / peakEnergy[ireg][ie][ires]);
1029  }
1030  ix++;
1031  }
1032  }
1033 
1034  //..Write the summary histograms to disk
1035  histfile->cd();
1036  peakSummary->Write();
1037  resolutionSummary->Write();
1038 
1039 
1040  //====================================================================================
1041  //====================================================================================
1042  //..Step 6. Finish up.
1043 
1044  //-----------------------------------------------------------------------------------
1045  //..Close histogram file
1046  histfile->Close();
1047 
1048  //------------------------------------------------------------------------
1049  //..Create and store payload
1050  ECLLeakageCorrections* leakagePayload = new ECLLeakageCorrections();
1051  leakagePayload->setlogEnergiesFwd(forwardVector);
1052  leakagePayload->setlogEnergiesBrl(barrelVector);
1053  leakagePayload->setlogEnergiesBwd(backwardVector);
1054  leakagePayload->setThetaCorrections(thetaCorrection);
1055  leakagePayload->setPhiCorrections(phiCorrection);
1056  leakagePayload->setnCrystalCorrections(nCrystalCorrection);
1057  saveCalibration(leakagePayload, "ECLLeakageCorrections");
1058 
1059  B2INFO("eclLeakageAlgorithm: successfully stored payload ECLLeakageCorrections");
1060  return c_OK;
1061 }
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_Failure
Failed =3 in Python.
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Singleton class to cache database objects.
Definition: DBStore.h:31
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:54
void setInitializeActive(bool active)
Setter for m_initializeActive.
Definition: DataStore.cc:94
DB object to store leakage corrections, including nCrys dependence
void setThetaCorrections(const TH2F &thetaCorrections)
Set the 2D histogram containing the theta corrections for each thetaID and energy.
void setlogEnergiesBrl(const std::vector< float > &logEnergiesBrl)
Set the vector of energies used to evaluate the leakage corrections in the barrel.
void setnCrystalCorrections(const TH2F &nCrystalCorrections)
Set the 2D histogram containing the nCrys corrections for each thetaID and energy.
void setlogEnergiesFwd(const std::vector< float > &logEnergiesFwd)
Set the vector of energies used to evaluate the leakage corrections in the forward endcap.
void setPhiCorrections(const TH2F &phiCorrections)
Set the 2D histogram containing the phi corrections for each thetaID and energy.
void setlogEnergiesBwd(const std::vector< float > &logEnergiesBwd)
Set the vector of energies used to evaluate the leakage corrections in the backward endcap.
float t_energyFrac
measured energy after nOptimal bias and peak corrections, divided by generated
float t_locationError
reconstructed minus generated position (cm)
int t_energyBin
generated energy point
int t_phiBin
binned location in phi relative to crystal edge
int t_phiMech
mechanical structure next to lower phi (0), upper phi (1), or neither (2)
int t_region
region of photon 0=forward 1=barrel 2=backward
int t_thetaBin
binned location in theta relative to crystal edge
int t_nCrys
number of crystals used to calculate energy
virtual EResult calibrate() override
Run algorithm.
double m_lowEnergyThreshold
Parameters to control fit procedure.
float t_origEnergyFrac
corrected energy at time of generation, divided by generated
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
bool construct(Args &&... params)
Construct an object of type T in this StoreObjPtr, using the provided constructor arguments.
Definition: StoreObjPtr.h:119
static DBStore & Instance()
Instance of a singleton DBStore.
Definition: DBStore.cc:28
void update()
Updates all objects that are outside their interval of validity.
Definition: DBStore.cc:79
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.