Belle II Software  release-06-02-00
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/dbobjects/ECLLeakageCorrections.h>
10 #include <framework/datastore/StoreObjPtr.h>
11 #include <framework/dataobjects/EventMetaData.h>
12 #include <framework/datastore/DataStore.h>
13 
14 
15 #include "TH1D.h"
16 #include "TF1.h"
17 #include "TTree.h"
18 #include "TFile.h"
19 #include "TDirectory.h"
20 #include "TProfile.h"
21 #include <iostream>
22 
23 using namespace std;
24 using namespace Belle2;
25 using namespace ECL;
26 using namespace Calibration;
27 
28 /**************************************************************************
29  * eclLeakageAlgorithm analyzes single photon MC to find the eclLeakage payload
30  * After initial setup, the position dependent corrections are found using
31  * photons with good reconstruction. This correction is then applied prior
32  * to finding the correction that depends on the number of crystals.
33  *
34  * Major steps in the process:
35  * 1. Set up (read in parameters and tree, define histogram binning)
36  * 2. Loop through data to histogram difference between true and reconstructed
37  * location. Determine cut for each thetaID/energy, used to select photons for
38  * all future steps.
39  * 3. Fill and fit histograms of normalized reconstructed energy, one per
40  * thetaID/energy. Novosibirsk fit parameter eta is floated in this fit, then
41  * fixed for fits of data with finer location binning. Peak (i.e. overall correction
42  * for this thetaID/energy) is used to normalize position-dependent correction.
43  * 4. Fill and fit histograms of normalized energy for each location (nominally 29
44  * locations in theta and phi per thetaID) to get the position-dependent correction.
45  * Nominally 69 x 8 x 29 x 3 = 48,024 histograms. 3 = theta, phi next to mech, phi not
46  * next to mech.
47  * 5. Fill and fit histograms of normalized reconstructed energy corrected for position
48  * dependent leakage, one per thetaID/energy. Resulting eta used for nCrys fits.
49  * 6. Fill and fit histograms of normalized energy corrected for location for each value
50  * of nCrys, for each thetaID / energy (nominally 69 x 8 x 21 = 11,592 histograms).
51  * Fix up corrections for those values of nCrys without sufficient statistics.
52  * 7. Pack payload quantities into vectors and histograms
53  * 8. Fill and fit summary and resolution histograms
54  * 9. Finish; close histogram file and and write payload
55 
56  **************************************************************************/
57 
58 
60 //..Function to get starting parameters for the Novo fit
61 std::vector<double> eclLeakageFitParameters(TH1F* h, const double& target)
62 {
63 
64  //..Find the smallest range of bins that contain at least "target" entries
65  double maxIntegral = h->Integral();
66  const int nBins = h->GetNbinsX();
67  int minLo = 1;
68  int minHi = nBins;
69 
70  //..Store a vector of integrals over histogram ranges so that I only need to
71  // call h->Integral nBins times, instead of nBins^2.
72  std::vector<double> intVector; // intVector[n] = sum of histogram bins [1,n+1] inclusive
73  intVector.push_back(h->GetBinContent(1));
74  for (int iLo = 2; iLo <= nBins; iLo++) {
75  double nextIntegral = intVector[iLo - 2] + h->GetBinContent(iLo);
76  intVector.push_back(nextIntegral);
77  }
78 
79  //..Now just look through all possible bin ranges to find the best one
80  for (int iLo = 2; iLo <= nBins; iLo++) {
81  for (int iHi = iLo; iHi <= nBins; iHi++) {
82 
83  //..sum[iLo, iHi] = sum[1, iHi] - sum[1,iLo-1] = intVector[iHi-1] - intVector[iLo-2]
84  double integral = intVector[iHi - 1] - intVector[iLo - 2];
85 
86  //..same number of bins, pick the pair with more entries
87  if ((integral > target and (iHi - iLo) < (minHi - minLo)) or
88  (integral > target and (iHi - iLo) == (minHi - minLo) and integral > maxIntegral)
89  ) {
90  minLo = iLo;
91  minHi = iHi;
92  maxIntegral = integral;
93  }
94  }
95  }
96 
97  //..Fit parameters are derived (in part) from this range of bins
98  const double lowE = h->GetBinLowEdge(minLo);
99  const double highE = h->GetBinLowEdge(minHi + 1);
100  const double peak = 0.5 * (lowE + highE);
101  const double sigma = 0.4 * (highE - lowE);
102  const double eta = 0.4; // typical value
103  const int nfitBins = (1 + minHi - minLo);
104 
105  std::vector<double> parameters;
106  parameters.push_back(peak);
107  parameters.push_back(sigma);
108  parameters.push_back(eta);
109  parameters.push_back(lowE);
110  parameters.push_back(highE);
111  parameters.push_back(nfitBins);
112 
113  return parameters;
114 }
115 
117 //..Novosibirsk; H. Ikeda et al. / NIM A 441 (2000) 401-426
118 double eclLeakageNovo(Double_t* x, Double_t* par)
119 {
120 
121  Double_t peak = par[1];
122  Double_t width = par[2];
123  Double_t sln4 = sqrt(log(4));
124  Double_t y = par[3] * sln4;
125  Double_t tail = -log(y + sqrt(1 + y * y)) / sln4;
126  Double_t qc = 0.;
127 
128  if (TMath::Abs(tail) < 1.e-7)
129  qc = 0.5 * TMath::Power(((x[0] - peak) / width), 2);
130  else {
131  double qa = tail * sqrt(log(4.));
132  double qb = sinh(qa) / qa;
133  double qx = (x[0] - peak) / width * qb;
134  double qy = 1. + tail * qx;
135 
136  if (qy > 1.E-7)
137  qc = 0.5 * (TMath::Power((log(qy) / tail), 2) + tail * tail);
138  else
139  qc = 15.0;
140  }
141 
142  return par[0] * exp(-qc);
143 }
144 
145 
147 //..Function to check for bad fits
148 // 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb, 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit
149 int eclLeakageFitQuality(const double& fitLo, const double& fitHi, const double& fitPeak, const double& fitSigma,
150  const double& fitdEta, const double& fitProb)
151 {
152  const double tolerance = 0.02;
153  const double redoFitProb = 1e-5;
154  const double minFitProb = 1e-9;
156  int fitStatus = 0;
157  if (fitProb < redoFitProb) {fitStatus = 1;}
158  if (fitProb < minFitProb) {fitStatus = 3;}
159  double tol = tolerance * (fitHi - fitLo);
160  if (fitPeak < (fitLo + tol) or fitPeak > (fitHi - tol)) {fitStatus = 4;}
161  if (fitSigma<tol or fitSigma>(fitHi - fitLo - tol)) {fitStatus = 5;}
162  double tolEta = tolerance;
163  if (fitdEta < tolEta) {fitStatus = 6;} // just use the tolerance for eta (note that it is eta - limit that is passed)
164  return fitStatus;
165 }
166 
168 eclLeakageAlgorithm::eclLeakageAlgorithm(): CalibrationAlgorithm("eclLeakageCollector")
169 {
171  "Generate payload ECLLeakageCorrection from single photon MC recorded by eclLeakageCollectorModule"
172  );
173 }
174 
175 
176 
179 {
180 
181  //====================================================================================
182  //====================================================================================
183  //..Step 1. Set up prior to first loop through data
184  B2INFO("starting eclLeakageAlgorithm");
185 
186  //-----------------------------------------------------------------------------------
187  //..Read in histograms created by the collector and fix normalization
188  auto inputParameters = getObjectPtr<TH1F>("inputParameters");
189  int lastBin = inputParameters->GetNbinsX();
190  double scale = inputParameters->GetBinContent(lastBin); // number of times inputParameters was filled
191  for (int ib = 1; ib < lastBin; ib++) {
192  double param = inputParameters->GetBinContent(ib);
193  inputParameters->SetBinContent(ib, param / scale);
194  inputParameters->SetBinError(ib, 0.);
195  }
196 
197  //..And write to disk.
198  TFile* histfile = new TFile("eclLeakageAlgorithm.root", "recreate");
199  inputParameters->Write();
200 
201  //..Parameters
202  const int nThetaID = 69;
203  const int nLeakReg = 3;
204  const int firstBarrelID = 13;
205  const int lastBarrelID = 58;
206  const int firstUsefulThID = 3;
207  const int lastUsefulThID = 66;
209  const int nPositions = (int)(inputParameters->GetBinContent(1) + 0.001);
210  const int nEnergies = (int)(inputParameters->GetBinContent(2) + 0.001);
211  const int nbinX = nEnergies * nThetaID;
213  const double etaMin = -5.;
214  const double etaMax = 5.;
217  //..Energies
218  float generatedE[nLeakReg][nEnergies]; // 3 regions forward barrel backward
219  int bin = 2; // bin 1 = nPositions, bin 2 = nEnergies, bin 3 = first energy
220  for (int ireg = 0; ireg < nLeakReg; ireg++) {
221  B2INFO("Generated energies for ireg = " << ireg);
222  for (int ie = 0; ie < nEnergies; ie++) {
223  bin++;
224  generatedE[ireg][ie] = inputParameters->GetBinContent(bin);
225  B2INFO(" " << ie << " " << generatedE[ireg][ie] << " GeV");
226  }
227  }
228  B2INFO("Low energy threshold = " << m_lowEnergyThreshold);
229  B2INFO("No nCrys threshold = " << m_noNCrysThreshold);
230 
231  //..Energy per thetaID (in MeV, for use in titles etc)
232  int iEnergiesMeV[nEnergies][nThetaID];
233  for (int thID = 0; thID < nThetaID; thID++) {
234  int ireg = 0;
235  if (thID >= firstBarrelID and thID <= lastBarrelID) {
236  ireg = 1;
237  } else if (thID > lastBarrelID) {
238  ireg = 2;
239  }
240  for (int ie = 0; ie < nEnergies; ie++) {
241  iEnergiesMeV[ie][thID] = (int)(1000.*generatedE[ireg][ie] + 0.5);
242  }
243  }
244 
245  //-----------------------------------------------------------------------------------
246  //..Bins for eFrac histograms (eFrac = uncorrected reconstructed E / E true)
247  const double eFracLo = 0.4; // low edge of eFrac histograms
248  const double eFracHi = 1.5; // high edge of eFrac histograms
249  int nEfracBins[nEnergies][nThetaID];
250  for (int thID = 0; thID < nThetaID; thID++) {
251  B2DEBUG(25, "eFrac nBins for thetaID " << thID);
252  for (int ie = 0; ie < nEnergies; ie++) {
253 
254  //..ballpark resolution
255  double res_squared = 0.0001 + 0.064 / iEnergiesMeV[ie][thID];
256 
257  //..Convert this to an even integer to get number of bins
258  double binNumOver2 = 3. / sqrt(res_squared);
259  int tempNBin = (int)(binNumOver2 + 0.5);
260  nEfracBins[ie][thID] = 2 * tempNBin;
261  B2DEBUG(25, " ie = " << ie << " E = " << iEnergiesMeV[ie][thID] << " nBins = " << nEfracBins[ie][thID]);
262  }
263  }
264 
265  //-----------------------------------------------------------------------------------
266  //..Set up Novosibirsk fit function
267  TString statusString[7] = {"good", "refit", "lowStat", "lowProb", "peakAtLimit", "sigmaAtLimit", "etaAtLimit"};
268  const double fracEnt[2] = {0.683, 0.5}; // fit range includes 68% or 50% of entries
269  const double minEntries = 100.; // don't use fits with fewer entries
270  const double minMaxBin = 50.; // rebin if max bin is below this value
271  const double highMaxBin = 300.; // can float eta if max bin is above this value
272 
273  TF1* func = new TF1("eclLeakageNovo", eclLeakageNovo, eFracLo, eFracHi, 4);
274  func->SetParNames("normalization", "peak", "effSigma", "eta");
275  func->SetLineColor(kRed);
276 
277  //-----------------------------------------------------------------------------------
278  //..Specify the TTree
279  auto tree = getObjectPtr<TTree>("tree");
280  tree->SetBranchAddress("cellID", &t_cellID);
281  tree->SetBranchAddress("thetaID", &t_thetaID);
282  tree->SetBranchAddress("region", &t_region);
283  tree->SetBranchAddress("thetaBin", &t_thetaBin);
284  tree->SetBranchAddress("phiBin", &t_phiBin);
285  tree->SetBranchAddress("phiMech", &t_phiMech);
286  tree->SetBranchAddress("energyBin", &t_energyBin);
287  tree->SetBranchAddress("nCrys", &t_nCrys);
288  tree->SetBranchAddress("energyFrac", &t_energyFrac);
289  tree->SetBranchAddress("origEnergyFrac", &t_origEnergyFrac);
290  tree->SetBranchAddress("locationError", &t_locationError);
291 
292  const int treeEntries = tree->GetEntries();
293  B2INFO("eclLeakageAlgorithm entries = " << treeEntries);
294 
295  //-----------------------------------------------------------------------------------
296  //..Get current payload to help validate the new payload
297 
298  //..Set event, run, exp number
299  const auto exprun = getRunList();
300  const int iEvt = 1;
301  const int iRun = exprun[0].second;
302  const int iExp = exprun[0].first;
305  evtPtr.registerInDataStore();
307  evtPtr.construct(iEvt, iRun, iExp);
308  DBStore& dbstore = DBStore::Instance();
309  dbstore.update();
310 
311  //..Existing payload
312  DBObjPtr<ECLLeakageCorrections> existingCorrections("ECLLeakageCorrections");
313  TH2F existingThetaCorrection = existingCorrections->getThetaCorrections();
314  existingThetaCorrection.SetName("existingThetaCorrection");
315  TH2F existingPhiCorrection = existingCorrections->getPhiCorrections();
316  existingPhiCorrection.SetName("existingPhiCorrection");
317  TH2F existingCrysCorrection = existingCorrections->getnCrystalCorrections();
318  existingCrysCorrection.SetName("existingnCrystalCorrection");
319 
320  //..Write out the correction histograms
321  histfile->cd();
322  existingThetaCorrection.Write();
323  existingPhiCorrection.Write();
324  existingCrysCorrection.Write();
325 
326  //====================================================================================
327  //====================================================================================
328  //..Step 2. First loop. Derive location cut for each energy and thetaID
329 
330  //-----------------------------------------------------------------------------------
331  //..Histograms of location error for each energy and thetaID
332  TH1F* locError[nEnergies][nThetaID];
333  TString name;
334  TString title;
335  for (int thID = 0; thID < nThetaID; thID++) {
336  for (int ie = 0; ie < nEnergies; ie++) {
337  name = "locError_" + std::to_string(ie) + "_" + std::to_string(thID);
338  title = "Location error " + to_string(iEnergiesMeV[ie][thID]) + " MeV thetaID " + std::to_string(thID) + ";location error (cm)";
339  locError[ie][thID] = new TH1F(name, title, 300, 0, 30.);
340  }
341  }
342 
343  //-----------------------------------------------------------------------------------
344  //..Loop over tree and histogram the location errors for useful thetaIDs
345  for (int i = 0; i < treeEntries; i++) {
346  tree->GetEntry(i);
347  if (t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0) {
348  locError[t_energyBin][t_thetaID]->Fill(t_locationError);
349  }
350  }
351 
352  //-----------------------------------------------------------------------------------
353  //..Cut is location where the distribution drops below 2.5% of peak after pedestal subtraction
354  const double minLocEntries = 49.5; // Write to disk if sufficient entries
355  const double peakFrac = 0.025; // look for distribution to drop to this fraction of peak value
356  const double startOfPed = 10.0001; // cm
357 
358  float maxLocCut[nEnergies][nThetaID]; // location cut for each energy and thetaID
359 
360  //..Summary histograms of cut values and fraction of events passing cut
361  TH1F* locCutSummary = new TH1F("locCutSummary", "location cut for each thetaID/energy; xBin = thetaID + ie*nTheta", nbinX, 0,
362  nbinX);
363  TH1F* locCutEfficiency = new TH1F("locCutEfficiency",
364  "Fraction of events passing location cut for each thetaID/energy; xBin = thetaID + ie*nTheta", nbinX, 0, nbinX);
365 
366  for (int thID = 0; thID < nThetaID; thID++) {
367  for (int ie = 0; ie < nEnergies; ie++) {
368  maxLocCut[ie][thID] = 0.;
369  if (locError[ie][thID]->Integral() < minLocEntries) {continue;}
370 
371  //..Pedestal per bin is average from 10 cm onwards
372  const int i10cm = locError[ie][thID]->GetXaxis()->FindBin(startOfPed);
373  const int ifinal = locError[ie][thID]->GetNbinsX();
374  double pedestal = locError[ie][thID]->Integral(i10cm, ifinal) / (ifinal + 1 - i10cm);
375  double threshold = pedestal + peakFrac * (locError[ie][thID]->GetMaximum() - pedestal);
376  const int ipeak = locError[ie][thID]->GetMaximumBin();
377 
378  //..Lower edge of 2nd consecutive bin below threshold
379  const int ix = 1 + thID + nThetaID * ie; // summary hist bin; +1 because first hist bin is 1, not 0
380  for (int ib = ipeak; ib <= ifinal; ib++) {
381  const double ibEnt = locError[ie][thID]->GetBinContent(ib);
382  const double ibMinus1 = locError[ie][thID]->GetBinContent(ib - 1);
383 
384  //..This is the cut location. Record in array and in summary histogram
385  if (ibEnt < threshold and ibMinus1 < threshold) {
386  maxLocCut[ie][thID] = locError[ie][thID]->GetBinLowEdge(ib);
387  locCutSummary->SetBinContent(ix, maxLocCut[ie][thID]);
388  locCutSummary->SetBinError(ix, 0);
389 
390  //..Efficiency, sort of
391  const double entries = locError[ie][thID]->GetEntries();
392  const double pass = locError[ie][thID]->Integral(1, ib - 1); // bin ib is not included
393  locCutEfficiency->SetBinContent(ix, pass / entries);
394  locCutEfficiency->SetBinError(ix, 0);
395 
396  break;
397  }
398  }
399  B2DEBUG(25, " thID " << thID << " E " << iEnergiesMeV[ie][thID] << " ped " << pedestal << " threshold " << threshold <<
400  " location cut " << maxLocCut[ie][thID]);
401 
402  //..Write this one to disk, after supplementing the title with cut and efficiency
403  title = locError[ie][thID]->GetTitle();
404  TString sCut;
405  const double locEff = locCutEfficiency->GetBinContent(ix);
406  sCut.Form(" cut %0.1f cm, eff %0.2f", maxLocCut[ie][thID], locEff);
407  title += sCut;
408  locError[ie][thID]->SetTitle(title);
409  histfile->cd();
410  locError[ie][thID]->Write();
411  }
412  }
413 
414  //..Write out the summary histogram
415  histfile->cd();
416  locCutSummary->Write();
417  locCutEfficiency->Write();
418 
419  //====================================================================================
420  //====================================================================================
421  //..Step 3. Second loop, fill histograms of e/eTrue for each energy and thetaID
422 
423  //-----------------------------------------------------------------------------------
424  //..One histogram of e/eTrue per energy per thetaID.
425  // Used to fix eta in subsequent fits, and to get overall correction for that thetaID
426  TH1F* hELabUncorr[nEnergies][nThetaID];
427  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
428  TString sthID = std::to_string(thID);
429  for (int ie = 0; ie < nEnergies; ie++) {
430  name = "hELabUncorr_" + std::to_string(ie) + "_" + sthID;
431  title = "eFrac " + to_string(iEnergiesMeV[ie][thID]) + " MeV thetaID " + sthID + ";E/Etrue";
432 
433  //..High statistics for these plots; use more bins
434  hELabUncorr[ie][thID] = new TH1F(name, title, 2 * nEfracBins[ie][thID], eFracLo, eFracHi);
435  }
436  }
437 
438  //-----------------------------------------------------------------------------------
439  //..Loop over events and store eFrac
440  for (int i = 0; i < treeEntries; i++) {
441  tree->GetEntry(i);
442 
443  //..Only events with good reconstruction
444  bool goodReco = t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0
445  and t_locationError < maxLocCut[t_energyBin][t_thetaID];
446  if (not goodReco) {continue;}
447 
448  //..Fill histogram for full thetaID
449  hELabUncorr[t_energyBin][t_thetaID]->Fill(t_energyFrac);
450  }
451 
452  //-----------------------------------------------------------------------------------
453  //..Fit each thetaID/energy histogram to get peak (overall correction) and eta (fixed
454  // in subsequent fits to individual locations).
455  float peakUncorr[nEnergies][nThetaID]; // store peak from each fit
456  float etaUncorr[nEnergies][nThetaID]; // store eta from each fit
457  std::vector<TString> failedELabUncorr; // names of hists with failed fits
458  std::vector<int> statusELabUncorr; // status of failed fits
459  int payloadStatus = 0; // Overall status of payload determination
460 
461  //..Record fit status
462  TH1F* statusOfhELabUncorr = new TH1F("statusOfhELabUncorr",
463  "status of hELabUncorr fits for each thetaID/energy. 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb, 4 = peakAtLimit, 5 = sigmaAtLimit",
464  6, 0,
465  6);
466 
467 
468  //..Fit each histogram
469  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
470  TString sthID = std::to_string(thID);
471  for (int ie = 0; ie < nEnergies; ie++) {
472  TH1F* hEnergy = (TH1F*)hELabUncorr[ie][thID];
473  double peak = -1.;
474  double eta = 0.;
475  int fitStatus = 2;
476  double entries = hEnergy->Integral();
477  int nIter = 0; // keep track of attempts to fit this histogram
478  double genE = iEnergiesMeV[ie][thID] / 1000.;
479  bool fitHist = entries > minEntries;
480 
481  //..Possibly iterate fit starting from this point
482  while (fitHist) {
483 
484  //..Fit parameters
485  double norm = hEnergy->GetMaximum();
486  double target = fracEnt[nIter] * entries; // fit range contains 68% or 50%
487  std::vector<double> startingParameters;// peak, sigma, eta, fitLow, fitHigh, nbins
488  startingParameters = eclLeakageFitParameters(hEnergy, target);
489  func->SetParameters(norm, startingParameters[0], startingParameters[1], startingParameters[2]);
490  func->SetParLimits(1, startingParameters[3], startingParameters[4]);
491  func->SetParLimits(2, 0., startingParameters[4] - startingParameters[3]);
492  func->SetParLimits(3, etaMin, etaMax);
493 
494  //..Fit
495  name = hEnergy->GetName();
496  B2DEBUG(40, "Fitting " << name.Data());
497  hEnergy->Fit(func, "LIEQ", "", startingParameters[3], startingParameters[4]);
498  peak = func->GetParameter(1);
499  double effSigma = func->GetParameter(2);
500  eta = func->GetParameter(3);
501  double prob = func->GetProb();
502 
503  //..Check fit quality 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb,
504  // 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit
505  double dEta = min((etaMax - eta), (eta - etaMin));
506  fitStatus = eclLeakageFitQuality(startingParameters[3], startingParameters[4], peak, effSigma, dEta, prob);
507 
508  //..If the fit probability is low, refit using a smaller range (fracEnt)
509  if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
510  nIter++;
511  } else {
512  fitHist = false;
513  }
514  }
515 
516  //-----------------------------------------------------------------------------------
517  //..Record failures and fit results.
518  // Mark payload as failed if energy is above low-energy threshold.
519  statusOfhELabUncorr->Fill(fitStatus + 0.000001);
520  if (fitStatus == 2 or fitStatus >= 4) {
521  statusELabUncorr.push_back(fitStatus);
522  failedELabUncorr.push_back(hEnergy->GetName());
523  if (genE > m_lowEnergyThreshold) {
524  payloadStatus = 1; // algorithm failed
525  } else {
526  peak = -1.; // fix up later
527  }
528  }
529  peakUncorr[ie][thID] = peak;
530  etaUncorr[ie][thID] = eta;
531 
532  //..Write to disk
533  if (entries > minEntries) {
534  histfile->cd();
535  hELabUncorr[ie][thID]->Write();
536  }
537  }
538  }
539 
540  //..Write out summary of fit status
541  histfile->cd();
542  statusOfhELabUncorr->Write();
543 
544  //-----------------------------------------------------------------------------------
545  //..Quit now if one of the high-energy fits failed, since we will not get a payload.
546  int nbadFit = statusELabUncorr.size();
547  if (nbadFit > 0) {B2ERROR("hELabUncorr fit failures (one histogram per energy/thetaID):");}
548  for (int ibad = 0; ibad < nbadFit; ibad++) {
549  int badStat = statusELabUncorr[ibad];
550  B2ERROR(" histogram " << failedELabUncorr[ibad].Data() << " status " << badStat << " " << statusString[badStat].Data());
551  }
552  if (payloadStatus != 0) {
553  B2ERROR("ecLeakageAlgorithm: fit to hELabUncorr failed. ");
554  histfile->Close();
555  return c_Failure;
556  }
557 
558  //-----------------------------------------------------------------------------------
559  //..Fix up any failed (low-energy) fits by using a neighbouring thetaID
560  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
561  for (int ie = 0; ie < nEnergies; ie++) {
562  if (peakUncorr[ie][thID] < 0.) {
563  if (thID > 40) {
564  for (int nextID = thID - 1; thID >= firstUsefulThID; thID--) {
565  if (peakUncorr[ie][nextID] > 0.) {
566  peakUncorr[ie][thID] = peakUncorr[ie][nextID];
567  break;
568  }
569  }
570  } else {
571  for (int nextID = thID + 1; thID <= lastUsefulThID; thID++) {
572  if (peakUncorr[ie][nextID] > 0.) {
573  peakUncorr[ie][thID] = peakUncorr[ie][nextID];
574  break;
575  }
576  }
577  }
578  }
579 
580  //..If we were unable to get a successful fit from a neighbour, give up
581  if (peakUncorr[ie][thID] < 0.) {
582  B2ERROR("ecLeakageAlgorithm: unable to correct hELabUncorr for thetaID " << thID << " ie " << ie);
583  histfile->Close();
584  return c_Failure;
585  }
586  }
587  }
588 
589  //====================================================================================
590  //====================================================================================
591  //..Step 4. Third loop, fill histograms of e/eTrue as a function of position.
592  // 29 locations in theta, 29 in phi.
593  // Crystals next to mechanical structure in phi are treated separately from
594  // crystals without mechanical structure.
595 
596  //-----------------------------------------------------------------------------------
597  //..Histograms to store the energy
598  const int nDir = 3;
599  TString dirName[nDir] = {"theta", "phiMech", "phiNoMech"};
600 
601  TH1F* eFracPosition[nEnergies][nThetaID][nDir][nPositions]; // the histograms
602  std::vector<TString> failedeFracPosition; // names of hists with failed fits
603  std::vector<int> statuseFracPosition; // status of failed fits
604 
605  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
606  TString sthID = std::to_string(thID);
607  for (int ie = 0; ie < nEnergies; ie++) {
608  TString sie = std::to_string(ie);
609  for (int idir = 0; idir < nDir; idir++) {
610  TString sidir = std::to_string(idir);
611  for (int ipos = 0; ipos < nPositions; ipos++) {
612  TString sipos = std::to_string(ipos);
613  name = "eFracPosition_" + sie + "_" + sthID + "_" + sidir + "_" + sipos;
614  title = "eFrac " + to_string(iEnergiesMeV[ie][thID]) + " MeV thetaID " + sthID + " " + dirName[idir] + " pos " + sipos +
615  "; E/Etrue";
616  eFracPosition[ie][thID][idir][ipos] = new TH1F(name, title, nEfracBins[ie][thID], eFracLo, eFracHi);
617  }
618  }
619  }
620  }
621 
622  //..And some summary histograms
623  TH1F* statusOfeFracPosition = new TH1F("statusOfeFracPosition",
624  "eFrac fit status: -2 noData, -1 lowE, 0 good, 1 redo fit, 2 lowStat, 3 lowProb, 4 peakAtLimit, 5 sigmaAtLimit", 8, -2,
625  6);
626  TH1F* probOfeFracPosition = new TH1F("probOfeFracPosition", "fit probability of eFrac fits for each position;probability", 100, 0,
627  1);
628  TH1F* maxOfeFracPosition = new TH1F("maxOfeFracPosition", "max entries of eFrac histograms;maximum bin content", 100, 0, 1000);
629 
630  //-----------------------------------------------------------------------------------
631  //..Loop over events and store eFrac
632  for (int i = 0; i < treeEntries; i++) {
633  tree->GetEntry(i);
634 
635  //..Only events with good reconstruction
636  bool goodReco = t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0
637  and t_locationError < maxLocCut[t_energyBin][t_thetaID];
638  if (not goodReco) {continue;}
639 
640  //..Theta location (idir = 0)
641  int idir = 0;
642  eFracPosition[t_energyBin][t_thetaID][idir][t_thetaBin]->Fill(t_energyFrac);
643 
644  //..phi location. Note that t_phiBin starts at mechanical edge, if there is one.
645  idir = t_phiMech + 1;
646  eFracPosition[t_energyBin][t_thetaID][idir][t_phiBin]->Fill(t_energyFrac);
647  }
648 
649  //-----------------------------------------------------------------------------------
650  //..Now fit many many histograms to get the position-dependent corrections
651  float positionCorrection[nEnergies][nThetaID][nDir][nPositions];
652  float positionCorrectionUnc[nEnergies][nThetaID][nDir][nPositions];
653  int nHistToFit = 0;
654 
655 
656  //..Temp histogram of position corrections
657  TH1F* thetaCorSummary = new TH1F("thetaCorSummary", "Theta dependent corrections;theta dependent correction", 100, 0.4, 1.4);
658  TH1F* phiCorSummary = new TH1F("phiCorSummary", "Phi dependent corrections;phi dependent correction", 100, 0.4, 1.4);
659 
660  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
661  for (int ie = 0; ie < nEnergies; ie++) {
662  double genE = iEnergiesMeV[ie][thID] / 1000.;
663  for (int idir = 0; idir < nDir; idir++) {
664  for (int ipos = 0; ipos < nPositions; ipos++) {
665  TH1F* hEnergy = (TH1F*)eFracPosition[ie][thID][idir][ipos];
666  if (hEnergy->Integral() > minEntries) {nHistToFit++;}
667 
668  //..Default peak from fit to full thetaID/energy = peakUncorr;
669  // correction = peak / sqrt(peakUncorr) = sqrt(peakUncorr)
670  double correction = sqrt(peakUncorr[ie][thID]);
671  double corrUnc = 0.05; // arbitrary uncertainty
672  double prob = -1.;
673  int fitStatus = 2; // low stats
674  if (genE < m_lowEnergyThreshold) {fitStatus = -1;} // low energy, don't fit, use default peak value
675  if (hEnergy->GetEntries() < 0.5) {fitStatus = -2;} // unused, eg barrel pos=2
676  int nIter = 0; // keep track of attempts to fit this histogram
677  double entries = hEnergy->Integral();
678  bool fitHist = entries > minEntries and genE >= m_lowEnergyThreshold;
679 
680  //..Possibly iterate fit starting from this point
681  while (fitHist) {
682  fitHist = false;
683 
684  //..Fit parameters
685  double target = fracEnt[nIter] * entries; // fit range contains 68%
686  std::vector<double> startingParameters;// peak, sigma, eta, fitLow, fitHigh, nbins
687  bool findParm = true;
688  while (findParm) {
689  startingParameters = eclLeakageFitParameters(hEnergy, target);
690 
691  //..Rebin if stats are low and we have enough DOF
692  if (hEnergy->GetMaximum()<minMaxBin and startingParameters[5]>11.9) {
693  hEnergy->Rebin(2);
694  } else {
695  findParm = false;
696  }
697  }
698  double norm = hEnergy->GetMaximum(); // max bin after rebinning
699  double fitLow = startingParameters[3];
700  double fitHigh = startingParameters[4];
701 
702  //..Eta from the fit to the full crystal
703  double etaFix = etaUncorr[ie][thID];
704  func->SetParameters(norm, startingParameters[0], startingParameters[1], etaFix);
705  func->SetParLimits(1, fitLow, fitHigh);
706  func->SetParLimits(2, 0., fitHigh - fitLow);
707 
708  //..Fix eta, unless really good statistics
709  if (hEnergy->GetMaximum() < highMaxBin) {
710  func->SetParLimits(3, etaFix, etaFix);
711  } else {
712  func->SetParLimits(3, etaMin, etaMax);
713  }
714 
715  //..Fit
716  name = hEnergy->GetName();
717  B2DEBUG(40, "Fitting " << name.Data());
718  hEnergy->Fit(func, "LIEQ", "", fitLow, fitHigh);
719  double peak = func->GetParameter(1);
720  double effSigma = func->GetParameter(2);
721  double eta = func->GetParameter(3);
722  prob = func->GetProb();
723 
724  //..Check fit quality 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb,
725  // 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit.
726  double dEta = min((etaMax - eta), (eta - etaMin));
727  fitStatus = eclLeakageFitQuality(fitLow, fitHigh, peak, effSigma, dEta, prob);
728 
729  //..If the fit probability is low, refit using a smaller range (fracEnt)
730  if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
731  nIter++;
732  fitHist = true;
733 
734  // Store correction except for peak or sigma at limit.
735  } else if (fitStatus <= 3) {
736 
737  //..Divide by sqrt(peak for full crystal) to get correct average
738  // when multiplying the theta and phi corrections
739  correction = peak / sqrt(peakUncorr[ie][thID]);
740  corrUnc = func->GetParError(1) / sqrt(peakUncorr[ie][thID]);
741  }
742  }
743 
744  //..Store the correction for this position
745  positionCorrection[ie][thID][idir][ipos] = correction;
746  positionCorrectionUnc[ie][thID][idir][ipos] = corrUnc;
747  statusOfeFracPosition->Fill(fitStatus + 0.00001);
748  probOfeFracPosition->Fill(prob);
749  maxOfeFracPosition->Fill(hEnergy->GetMaximum());
750 
751  //..Summary histogram for successful fits. Note that fitStatus<0 also has prob<0.
752  if (prob > 0 and fitStatus <= 3) {
753  if (idir == 0) {
754  thetaCorSummary->Fill(correction);
755  } else {
756  phiCorSummary->Fill(correction);
757  }
758  }
759 
760  //..Record the failed fit details. We may still use the correction.
761  if (fitStatus >= 2) {
762  statuseFracPosition.push_back(fitStatus);
763  failedeFracPosition.push_back(hEnergy->GetName());
764  histfile->cd();
765  hEnergy->Write();
766  }
767  }
768  }
769  }
770  }
771 
772  //-----------------------------------------------------------------------------------
773  //..Summarize position fit results
774  B2INFO("eclLeakageAlgorithm: " << nHistToFit << " eFracPosition histograms to fit");
775  nbadFit = statuseFracPosition.size();
776  B2INFO("eFracPosition failed fits: " << nbadFit);
777  for (int ibad = 0; ibad < nbadFit; ibad++) {
778  int badStat = statuseFracPosition[ibad];
779  B2ERROR(" histogram " << failedeFracPosition[ibad].Data() << " status " << badStat << " " << statusString[badStat].Data());
780  }
781 
782  //..Write to disk
783  histfile->cd();
784  statusOfeFracPosition->Write();
785  probOfeFracPosition->Write();
786  maxOfeFracPosition->Write();
787  thetaCorSummary->Write();
788  phiCorSummary->Write();
789 
790  //====================================================================================
791  //====================================================================================
792  //..Step 5. Fourth loop, fill and fit histograms of ePos for each thetaID/energy,
793  // where ePos = normalized energy after position correction.
794 
795  //-----------------------------------------------------------------------------------
796  //..One histogram of ePos per energy per thetaID.
797  TH1F* hEPos[nEnergies][nThetaID];
798  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
799  TString sthID = std::to_string(thID);
800  for (int ie = 0; ie < nEnergies; ie++) {
801  name = "hEPos_" + std::to_string(ie) + "_" + sthID;
802  title = "ePos " + to_string(iEnergiesMeV[ie][thID]) + " MeV thetaID " + sthID + ";Position corrected E/Etrue";
803 
804  //..High statistics for these plots; use more bins
805  hEPos[ie][thID] = new TH1F(name, title, 2 * nEfracBins[ie][thID], eFracLo, eFracHi);
806  }
807  }
808 
809  //-----------------------------------------------------------------------------------
810  //..Loop over events and store ePos
811  for (int i = 0; i < treeEntries; i++) {
812  tree->GetEntry(i);
813 
814  //..Only events with good reconstruction
815  bool goodReco = t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0
816  and t_locationError < maxLocCut[t_energyBin][t_thetaID];
817  if (not goodReco) {continue;}
818 
819  //..Position-dependent leakage corrections
820  int idir = 0;
821  float thetaPosCor = positionCorrection[t_energyBin][t_thetaID][idir][t_thetaBin];
822 
823  idir = t_phiMech + 1;
824  float phiPosCor = positionCorrection[t_energyBin][t_thetaID][idir][t_phiBin];
825 
826  float ePos = t_energyFrac / thetaPosCor / phiPosCor;
827 
828  //..Fill histogram for full thetaID
829  hEPos[t_energyBin][t_thetaID]->Fill(ePos);
830  }
831 
832  //-----------------------------------------------------------------------------------
833  //..Fit each thetaID/energy histogram to get eta (fixed for later nCrys fits)
834  float etaEpos[nEnergies][nThetaID]; // store eta from each fit
835  std::vector<TString> failedEPos; // names of hists with failed fits
836  std::vector<int> statusEPos; // status of failed fits
837 
838  //..Record fit status
839  TH1F* statusOfhEPos = new TH1F("statusOfhEPos",
840  "status of hEPos fits for each thetaID/energy: -2 noData, -1 lowE, 0 good, 1 redo fit, 2 lowStat, 3 lowProb, 4 peakAtLimit, 5 sigmaAtLimit",
841  8, -2,
842  6);
843 
844  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
845  TString sthID = std::to_string(thID);
846  for (int ie = 0; ie < nEnergies; ie++) {
847  double genE = iEnergiesMeV[ie][thID] / 1000.;
848  TH1F* hEnergy = (TH1F*)hEPos[ie][thID];
849  double eta = 0.4; // arbitrary but typical
850  int fitStatus = 2; // 2 = low stats
851  if (genE < m_lowEnergyThreshold or genE < m_noNCrysThreshold) {fitStatus = -1;} // low energy, don't fit
852  double entries = hEnergy->Integral();
853  int nIter = 0; // keep track of attempts to fit this histogram
854  bool fitHist = entries > minEntries and genE >= m_lowEnergyThreshold and genE >= m_noNCrysThreshold;
855 
856  //..Possibly iterate fit starting from this point
857  while (fitHist) {
858 
859  //..Fit parameters
860  double norm = hEnergy->GetMaximum();
861  double target = fracEnt[nIter] * entries; // fit range contains 68% or 50%
862  std::vector<double> startingParameters;// peak, sigma, eta, fitLow, fitHigh, nbins
863  startingParameters = eclLeakageFitParameters(hEnergy, target);
864  func->SetParameters(norm, startingParameters[0], startingParameters[1], startingParameters[2]);
865  func->SetParLimits(1, startingParameters[3], startingParameters[4]);
866  func->SetParLimits(2, 0., startingParameters[4] - startingParameters[3]);
867  func->SetParLimits(3, etaMin, etaMax);
868 
869  //..Fit
870  name = hEnergy->GetName();
871  B2DEBUG(40, "Fitting " << name.Data());
872  hEnergy->Fit(func, "LIEQ", "", startingParameters[3], startingParameters[4]);
873  double peak = func->GetParameter(1);
874  double effSigma = func->GetParameter(2);
875  eta = func->GetParameter(3);
876  double prob = func->GetProb();
877 
878  //..Check fit quality 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb,
879  // 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit
880  double dEta = min((etaMax - eta), (eta - etaMin));
881  fitStatus = eclLeakageFitQuality(startingParameters[3], startingParameters[4], peak, effSigma, dEta, prob);
882 
883  //..If the fit probability is low, refit using a smaller range (fracEnt)
884  if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
885  nIter++;
886  } else {
887  fitHist = false;
888  }
889  }
890 
891  //-----------------------------------------------------------------------------------
892  //..Record failures and fit results.
893  // Mark payload as failed if low stats or a fit at limit.
894  statusOfhEPos->Fill(fitStatus + 0.00001);
895  nbadFit = statusEPos.size();
896  if (nbadFit > 0) {B2ERROR("hEPos fit failures (one histogram per energy/thetaID after position correction):");}
897  if (fitStatus >= 2) {
898  statusEPos.push_back(fitStatus);
899  failedEPos.push_back(hEnergy->GetName());
900  if (fitStatus == 2 or fitStatus >= 4) { payloadStatus = 1; }
901  }
902  etaEpos[ie][thID] = eta; // value for low-energy cases is never used
903 
904  //..Write to disk
905  if (entries > minEntries) {
906  histfile->cd();
907  hEPos[ie][thID]->Write();
908  }
909  }
910  }
911 
912  //..Write status summary to disk
913  histfile->cd();
914  statusOfhEPos->Write();
915 
916  if (payloadStatus != 0) {
917  B2ERROR("ecLeakageAlgorithm: fit to hEPos failed. ");
918  histfile->Close();
919  return c_Failure;
920  }
921 
922 
923  //====================================================================================
924  //====================================================================================
925  //..Step 6. Fifth loop, fill and fit histograms of ePos for each nCrys, for each
926  // thetaID/energy, where ePos = normalized energy after position correction.
927  const int maxN = 21; // maximum number of crystals in a cluster
928 
929  //-----------------------------------------------------------------------------------
930  //..Histograms of ePos
931  TH1F* ePosnCry[nEnergies][nThetaID][maxN + 1]; // +1 to include 0
932  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
933  TString sthID = std::to_string(thID);
934  for (int ie = 0; ie < nEnergies; ie++) {
935  TString sie = std::to_string(ie);
936  for (int in = 0; in <= maxN; in++) {
937  name = "ePosnCry_" + sie + "_" + sthID + "_" + std::to_string(in);
938  title = "ePos " + to_string(iEnergiesMeV[ie][thID]) + " MeV thetaID " + sthID + " nCrys " + std::to_string(
939  in) + "; corrected E/Etrue";
940  ePosnCry[ie][thID][in] = new TH1F(name, title, nEfracBins[ie][thID], eFracLo, eFracHi);
941  }
942  }
943  }
944 
945  //-----------------------------------------------------------------------------------
946  //..Loop over events and store ePos
947  for (int i = 0; i < treeEntries; i++) {
948  tree->GetEntry(i);
949 
950  //..Only events with good reconstruction
951  bool goodReco = t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0
952  and t_locationError < maxLocCut[t_energyBin][t_thetaID];
953  if (not goodReco) {continue;}
954 
955  //..Position-dependent leakage corrections
956  int idir = 0;
957  float thetaPosCor = positionCorrection[t_energyBin][t_thetaID][idir][t_thetaBin];
958 
959  idir = t_phiMech + 1;
960  float phiPosCor = positionCorrection[t_energyBin][t_thetaID][idir][t_phiBin];
961 
962  float ePos = t_energyFrac / thetaPosCor / phiPosCor;
963  ePosnCry[t_energyBin][t_thetaID][t_nCrys]->Fill(ePos);
964  }
965 
966  //-----------------------------------------------------------------------------------
967  //..Fit ePos histograms for each value of nCrys. Many of these will not have enough stats
968  float nCrysCorrection[nEnergies][nThetaID][maxN + 1];
969  float nCrysCorrectionUnc[nEnergies][nThetaID][maxN + 1];
970 
971  //..Keep track of the peak nCrys for each thetaID/energy. Use this later to fix up
972  // values of nCrys without a successful fit
973  int maxNCry[nEnergies][nThetaID];
974 
975  //..Keep track of fit status
976  TH1F* statusOfePosnCry = new TH1F("statusOfePosnCry",
977  "ePosnCry fit status: -2 noData, -1 lowE, 0 good, 1 redo fit, 2 lowStat, 3 lowProb, 4 peakAtLimit, 5 sigmaAtLimit", 8, -2,
978  6);
979 
980  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
981  for (int ie = 0; ie < nEnergies; ie++) {
982  double genE = iEnergiesMeV[ie][thID] / 1000.;
983  double maxIntegral = 0;
984  maxNCry[ie][thID] = 0;
985  for (int in = 0; in <= maxN; in++) {
986 
987  TH1F* hEnergy = (TH1F*)ePosnCry[ie][thID][in];
988  if (hEnergy->Integral() > minEntries) {nHistToFit++;}
989 
990  //..Defaults for failed fits. Will fix up in the next step.
991  double correction = -1.;
992  double corrUnc = 0.05; // arbitrary uncertainty
993  int nIter = 0; // keep track of attempts to fit this histogram
994  double entries = hEnergy->Integral();
995  int fitStatus = -2; // low stats
996  if (genE < m_lowEnergyThreshold or genE < m_noNCrysThreshold) {
997  fitStatus = -1; // low energy
998  correction = 1.; // for low energy, nCrys correction = 1
999  corrUnc = 0.;
1000  }
1001  bool fitHist = entries > minEntries and genE >= m_lowEnergyThreshold and genE >= m_noNCrysThreshold;
1002 
1003  //..nCrys with most entries
1004  if (entries > maxIntegral) {
1005  maxIntegral = entries;
1006  maxNCry[ie][thID] = in;
1007  }
1008 
1009  //..Possibly iterate fit starting from this point
1010  while (fitHist) {
1011  fitHist = false;
1012 
1013  //..Fit parameters
1014  double target = fracEnt[nIter] * entries; // fit range contains 68%
1015  std::vector<double> startingParameters;// peak, sigma, eta, fitLow, fitHigh, nbins
1016  bool findParm = true;
1017  while (findParm) {
1018  startingParameters = eclLeakageFitParameters(hEnergy, target);
1019 
1020  //..Rebin if stats are low and we have enough DOF
1021  if (hEnergy->GetMaximum()<minMaxBin and startingParameters[5]>11.9) {
1022  hEnergy->Rebin(2);
1023  } else {
1024  findParm = false;
1025  }
1026  }
1027  double norm = hEnergy->GetMaximum(); // max bin after rebinning
1028  double fitLow = startingParameters[3];
1029  double fitHigh = startingParameters[4];
1030 
1031  //..Eta from the fit to the full crystal after position-dependent leakage correction
1032  double etaFix = etaEpos[ie][thID];
1033  func->SetParameters(norm, startingParameters[0], startingParameters[1], etaFix);
1034  func->SetParLimits(1, fitLow, fitHigh);
1035  func->SetParLimits(2, 0., fitHigh - fitLow);
1036 
1037  //..Fix eta, unless really good statistics
1038  if (hEnergy->GetMaximum() < highMaxBin) {
1039  func->SetParLimits(3, etaFix, etaFix);
1040  } else {
1041  func->SetParLimits(3, etaMin, etaMax);
1042  }
1043 
1044  //..Fit
1045  name = hEnergy->GetName();
1046  B2DEBUG(40, "Fitting " << name.Data());
1047  hEnergy->Fit(func, "LIEQ", "", fitLow, fitHigh);
1048  double peak = func->GetParameter(1);
1049  double effSigma = func->GetParameter(2);
1050  double eta = func->GetParameter(3);
1051  double prob = func->GetProb();
1052 
1053  //..Check fit quality 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb,
1054  // 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit.
1055  double dEta = min((etaMax - eta), (eta - etaMin));
1056  fitStatus = eclLeakageFitQuality(fitLow, fitHigh, peak, effSigma, dEta, prob);
1057 
1058  //..If the fit probability is low, refit using a smaller range (fracEnt)
1059  if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
1060  nIter++;
1061  fitHist = true;
1062 
1063  // Store correction except for peak or sigma at limit.
1064  } else if (fitStatus <= 3) {
1065  correction = peak;
1066  corrUnc = func->GetParError(1);
1067  }
1068  }
1069 
1070  //..Store the correction for this position
1071  statusOfePosnCry->Fill(fitStatus + 0.00001);
1072  nCrysCorrection[ie][thID][in] = correction;
1073  nCrysCorrectionUnc[ie][thID][in] = corrUnc;
1074  }
1075 
1076  //..Write histogram for maximum nCrys to disk
1077  histfile->cd();
1078  int highestN = maxNCry[ie][thID];
1079  ePosnCry[ie][thID][highestN]->Write();
1080  }
1081  }
1082 
1083  //..Write out fit status summary
1084  histfile->cd();
1085  statusOfePosnCry->Write();
1086 
1087  //-----------------------------------------------------------------------------------
1088  //..Now fix up corrections for nCrys values that did not have a successful fit.
1089  // The logic is that if "n" is not successful, then use the correction from either
1090  // n+1 or n-1, which ever is closer to the most likely value of nCrys
1091  for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
1092  for (int ie = 0; ie < nEnergies; ie++) {
1093 
1094  //..If the most likely value of nCrys does not have a successful correction,
1095  // mark the algorithm as failed. Recall that low-energy cases always have correction = 1.
1096  int highestN = maxNCry[ie][thID];
1097  if (nCrysCorrection[ie][thID][highestN] <= 0) {
1098  B2ERROR("eclLeakageAlgorithm: no nCrys correction for ie " << ie << " thetaID " << thID);
1099  histfile->Close();
1100  return c_Failure;
1101  }
1102 
1103  //..Start at the highest and work upwards
1104  for (int in = highestN + 1; in <= maxN; in++) {
1105  if (nCrysCorrection[ie][thID][in] < 0) {
1106  nCrysCorrection[ie][thID][in] = nCrysCorrection[ie][thID][in - 1];
1107  nCrysCorrectionUnc[ie][thID][in] = nCrysCorrectionUnc[ie][thID][in - 1];
1108  }
1109  }
1110 
1111  //..Start at the highest and move downwards
1112  for (int in = highestN - 1; in >= 0; in--) {
1113  if (nCrysCorrection[ie][thID][in] < 0) {
1114  nCrysCorrection[ie][thID][in] = nCrysCorrection[ie][thID][in + 1];
1115  nCrysCorrectionUnc[ie][thID][in] = nCrysCorrectionUnc[ie][thID][in + 1];
1116  }
1117  }
1118  }
1119  }
1120 
1121  //====================================================================================
1122  //====================================================================================
1123  //..Step 7. Store quantities needed for the payloads
1124 
1125  //-----------------------------------------------------------------------------------
1126  //..First, we need to fix up the thetaID that have been so far missed.
1127  // Just copy the values from the first or last thetaID for which corrections were found.
1128 
1129  //..ThetaID before the first useful one
1130  for (int thID = firstUsefulThID - 1; thID >= 0; thID--) {
1131  for (int ie = 0; ie < nEnergies; ie++) {
1132 
1133  //..Position dependent
1134  for (int idir = 0; idir < 3; idir++) {
1135  for (int ipos = 0; ipos < nPositions; ipos++) {
1136  positionCorrection[ie][thID][idir][ipos] = positionCorrection[ie][thID + 1][idir][ipos];
1137  positionCorrectionUnc[ie][thID][idir][ipos] = positionCorrectionUnc[ie][thID + 1][idir][ipos];
1138  }
1139  }
1140 
1141  //..nCrys
1142  for (int in = 0; in <= maxN; in++) {
1143  nCrysCorrection[ie][thID][in] = nCrysCorrection[ie][thID + 1][in];
1144  nCrysCorrectionUnc[ie][thID][in] = nCrysCorrectionUnc[ie][thID + 1][in];
1145  }
1146  }
1147  }
1148 
1149  //..ThetaID beyond last useful one
1150  for (int thID = lastUsefulThID + 1; thID < nThetaID; thID++) {
1151  for (int ie = 0; ie < nEnergies; ie++) {
1152 
1153  //..Position dependent
1154  for (int idir = 0; idir < 3; idir++) {
1155  for (int ipos = 0; ipos < nPositions; ipos++) {
1156  positionCorrection[ie][thID][idir][ipos] = positionCorrection[ie][thID - 1][idir][ipos];
1157  positionCorrectionUnc[ie][thID][idir][ipos] = positionCorrectionUnc[ie][thID - 1][idir][ipos];
1158  }
1159  }
1160 
1161  //..nCrys
1162  for (int in = 0; in <= maxN; in++) {
1163  nCrysCorrection[ie][thID][in] = nCrysCorrection[ie][thID - 1][in];
1164  nCrysCorrectionUnc[ie][thID][in] = nCrysCorrectionUnc[ie][thID - 1][in];
1165  }
1166  }
1167  }
1168 
1169  //-----------------------------------------------------------------------------------
1170  //..std::vectors of generated energies
1171  std::vector<float> forwardVector;
1172  std::vector<float> barrelVector;
1173  std::vector<float> backwardVector;
1174  for (int ie = 0; ie < nEnergies; ie++) {
1175  forwardVector.push_back(log(generatedE[0][ie]));
1176  barrelVector.push_back(log(generatedE[1][ie]));
1177  backwardVector.push_back(log(generatedE[2][ie]));
1178  }
1179 
1180  //..Store in array format for validation studies
1181  float leakLogE[3][12] = {};
1182  for (int ie = 0; ie < nEnergies; ie++) {
1183  leakLogE[0][ie] = forwardVector[ie];
1184  leakLogE[1][ie] = barrelVector[ie];
1185  leakLogE[2][ie] = backwardVector[ie];
1186  }
1187 
1188  //-----------------------------------------------------------------------------------
1189  //..Position dependent corrections
1190 
1191  //..2D histogram of theta-dependent corrections
1192  TH2F thetaCorrection("thetaCorrection", "Theta correction vs bin;bin = thetaID + 69*energyBin;local theta bin", nbinX, 0, nbinX,
1193  nPositions, 0, nPositions);
1194  int ix = 0;
1195  for (int ie = 0; ie < nEnergies; ie++) {
1196  for (int thID = 0; thID < nThetaID; thID++) {
1197  ix++;
1198  for (int ipos = 0; ipos < nPositions; ipos++) {
1199  int iy = ipos + 1;
1200  float correction = positionCorrection[ie][thID][0][ipos];
1201  float corrUnc = positionCorrectionUnc[ie][thID][0][ipos];
1202  thetaCorrection.SetBinContent(ix, iy, correction);
1203  thetaCorrection.SetBinError(ix, iy, corrUnc);
1204  }
1205  }
1206  }
1207 
1208  //..2D histogram of phi-dependent corrections. Twice as many x bins;
1209  // one set for crystals next to mechanical structure, one for otherwise.
1210  TH2F phiCorrection("phiCorrection", "Phi correction vs bin;bin = thetaID + 69*energyBin;local phi bin", 2 * nbinX, 0, 2 * nbinX,
1211  nPositions, 0, nPositions);
1212  ix = 0;
1213  for (int ie = 0; ie < nEnergies; ie++) {
1214  for (int thID = 0; thID < nThetaID; thID++) {
1215  ix++;
1216  for (int ipos = 0; ipos < nPositions; ipos++) {
1217  int iy = ipos + 1;
1218 
1219  //..First set of corrections
1220  float correction = positionCorrection[ie][thID][1][ipos];
1221  float corrUnc = positionCorrectionUnc[ie][thID][1][ipos];
1222  phiCorrection.SetBinContent(ix, iy, correction);
1223  phiCorrection.SetBinError(ix, iy, corrUnc);
1224 
1225  //..Second set
1226  correction = positionCorrection[ie][thID][2][ipos];
1227  corrUnc = positionCorrectionUnc[ie][thID][2][ipos];
1228  phiCorrection.SetBinContent(ix + nbinX, iy, correction);
1229  phiCorrection.SetBinError(ix + nbinX, iy, corrUnc);
1230  }
1231  }
1232  }
1233 
1234  //-----------------------------------------------------------------------------------
1235  //..nCrys dependent corrections
1236  TH2F nCrystalCorrection("nCrystalCorrection", "nCrys correction vs bin;bin = thetaID + 69*energyBin;nCrys", nbinX, 0, nbinX,
1237  maxN + 1, 0, maxN + 1);
1238 
1239  ix = 0;
1240  for (int ie = 0; ie < nEnergies; ie++) {
1241  for (int thID = 0; thID < nThetaID; thID++) {
1242  ix++;
1243  for (int in = 0; in <= maxN; in++) {
1244  int iy = in + 1;
1245  float correction = nCrysCorrection[ie][thID][in];
1246  float corrUnc = nCrysCorrectionUnc[ie][thID][in];
1247  nCrystalCorrection.SetBinContent(ix, iy, correction);
1248  nCrystalCorrection.SetBinError(ix, iy, corrUnc);
1249  }
1250  }
1251  }
1252 
1253 
1254  //====================================================================================
1255  //====================================================================================
1256  //..Step 8. Diagnostic histograms
1257 
1258  //..Start by storing the payload histograms
1259  histfile->cd();
1260  thetaCorrection.Write();
1261  phiCorrection.Write();
1262  nCrystalCorrection.Write();
1263 
1264  //-----------------------------------------------------------------------------------
1265  //..One histogram of new and original reconstructed energy after leakage correction
1266  // per generated energy per region. Also uncorrected.
1267  const int nResType = 5;
1268  TString resName[nResType] = {"Uncorrected", "Original", "Corrected no nCrys", "Corrected measured", "Corrected true"};
1269  TString regName[nLeakReg] = {"forward", "barrel", "backward"};
1270  TH1F* energyResolution[nLeakReg][nEnergies][nResType];
1271 
1272  //..Base number of bins on a typical thetaID for each region
1273  int thIDReg[nLeakReg];
1274  thIDReg[0] = 9;
1275  thIDReg[1] = 35;
1276  thIDReg[2] = 61;
1277 
1278  for (int ireg = 0; ireg < nLeakReg; ireg++) {
1279  for (int ie = 0; ie < nEnergies; ie++) {
1280 
1281  //..Titles and bins for this region and energy
1282  TString stireg = std::to_string(ireg);
1283  TString stie = std::to_string(ie);
1284  TString stieName = std::to_string(iEnergiesMeV[ie][thIDReg[ireg]]);
1285  int nbinReg = 20 * nEfracBins[ie][thIDReg[ireg]];
1286 
1287  for (int ires = 0; ires < nResType; ires++) {
1288  name = "energyResolution_" + stireg + "_" + stie + "_" + std::to_string(ires);
1289  title = resName[ires] + " energy, " + stieName + " MeV, " + regName[ireg] + ";Reconstructed E/Etrue";
1290  energyResolution[ireg][ie][ires] = new TH1F(name, title, nbinReg, eFracLo, eFracHi);
1291  }
1292  }
1293  }
1294 
1295  //-----------------------------------------------------------------------------------
1296  //..Loop over events and store leakage-corrected energies
1297  for (int i = 0; i < treeEntries; i++) {
1298  tree->GetEntry(i);
1299 
1300  //..Only events with good reconstruction
1301  bool goodReco = t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0
1302  and t_locationError < maxLocCut[t_energyBin][t_thetaID];
1303  if (not goodReco) {continue;}
1304 
1305  //-----------------------------------------------------------------------------------
1306  //..Corrections using true energy
1307 
1308  //..Position-dependent leakage corrections using true energy
1309  int idir = 0;
1310  float thetaPosCor = positionCorrection[t_energyBin][t_thetaID][idir][t_thetaBin];
1311 
1312  idir = t_phiMech + 1;
1313  float phiPosCor = positionCorrection[t_energyBin][t_thetaID][idir][t_phiBin];
1314  float posCor = thetaPosCor * phiPosCor;
1315 
1316  //..nCrystal-dependent leakage correction using true energy
1317  float nCrysCor = nCrysCorrection[t_energyBin][t_thetaID][t_nCrys];
1318  float corrTrue = posCor * nCrysCor;
1319 
1320  //-----------------------------------------------------------------------------------
1321  //..Find position-dependent correction (only) using measured energy. The correction is
1322  // a function of corrected energy, so will need to iterate
1323  float corrNonCrys = 0.96; // typical correction as starting point
1324  for (int iter = 0; iter < 2; iter++) {
1325 
1326 
1327  //..Energy points that bracket this value
1328  float energyRaw = t_energyFrac * generatedE[t_region][t_energyBin];
1329  float logEnergy = log(energyRaw / corrNonCrys);
1330  int ie0 = 0; // lower energy point
1331  if (logEnergy < leakLogE[t_region][0]) {
1332  ie0 = 0;
1333  } else if (logEnergy > leakLogE[t_region][nEnergies - 1]) {
1334  ie0 = nEnergies - 2;
1335  } else {
1336  while (logEnergy > leakLogE[t_region][ie0 + 1]) {ie0++;}
1337  }
1338 
1339  //..Corrections from lower and upper energy points
1340  float cor0 = positionCorrection[ie0][t_thetaID][0][t_thetaBin] * positionCorrection[ie0][t_thetaID][t_phiMech + 1][t_phiBin];
1341  float cor1 = positionCorrection[ie0 + 1][t_thetaID][0][t_thetaBin] * positionCorrection[ie0 + 1][t_thetaID][t_phiMech +
1342  1][t_phiBin];
1343 
1344  //..Interpolate in logE
1345  corrNonCrys = cor0 + (cor1 - cor0) * (logEnergy - leakLogE[t_region][ie0]) / (leakLogE[t_region][ie0 + 1] -
1346  leakLogE[t_region][ie0]);
1347  }
1348 
1349  //-----------------------------------------------------------------------------------
1350  //..Find correction using measured energy. The correction is a function of corrected
1351  // energy, so will need to iterate
1352  float corrMeasured = 0.96; // typical correction as starting point
1353  for (int iter = 0; iter < 2; iter++) {
1354 
1355 
1356  //..Energy points that bracket this value
1357  float energyRaw = t_energyFrac * generatedE[t_region][t_energyBin];
1358  float logEnergy = log(energyRaw / corrMeasured);
1359  int ie0 = 0; // lower energy point
1360  if (logEnergy < leakLogE[t_region][0]) {
1361  ie0 = 0;
1362  } else if (logEnergy > leakLogE[t_region][nEnergies - 1]) {
1363  ie0 = nEnergies - 2;
1364  } else {
1365  while (logEnergy > leakLogE[t_region][ie0 + 1]) {ie0++;}
1366  }
1367 
1368  //..Corrections from lower and upper energy points
1369  float cor0 = positionCorrection[ie0][t_thetaID][0][t_thetaBin] * positionCorrection[ie0][t_thetaID][t_phiMech + 1][t_phiBin] *
1370  nCrysCorrection[ie0][t_thetaID][t_nCrys];
1371  float cor1 = positionCorrection[ie0 + 1][t_thetaID][0][t_thetaBin] * positionCorrection[ie0 + 1][t_thetaID][t_phiMech +
1372  1][t_phiBin] * nCrysCorrection[ie0 + 1][t_thetaID][t_nCrys];
1373 
1374  //..Interpolate in logE
1375  corrMeasured = cor0 + (cor1 - cor0) * (logEnergy - leakLogE[t_region][ie0]) / (leakLogE[t_region][ie0 + 1] -
1376  leakLogE[t_region][ie0]);
1377  }
1378 
1379  //-----------------------------------------------------------------------------------
1380  //..Fill the histograms
1381  energyResolution[t_region][t_energyBin][0]->Fill(t_energyFrac); // uncorrected
1382  energyResolution[t_region][t_energyBin][1]->Fill(t_origEnergyFrac); // original payload
1383  energyResolution[t_region][t_energyBin][2]->Fill(t_energyFrac / corrNonCrys); // no nCrys, measured energy
1384  energyResolution[t_region][t_energyBin][3]->Fill(t_energyFrac / corrMeasured); // corrected, measured energy
1385  energyResolution[t_region][t_energyBin][4]->Fill(t_energyFrac / corrTrue); // corrected, true energy
1386  }
1387 
1388  //-----------------------------------------------------------------------------------
1389  //..Fit each histogram to find peak, and extract resolution.
1390 
1391  //..Store the peak and resolution values for each histogram
1392  float peakEnergy[nLeakReg][nEnergies][nResType];
1393  float energyRes[nLeakReg][nEnergies][nResType];
1394 
1395  //..Loop over the histograms to be fit
1396  for (int ireg = 0; ireg < nLeakReg; ireg++) {
1397  for (int ie = 0; ie < nEnergies; ie++) {
1398 
1399  //..Base the resolution on the number of entries in the corrected plot for all 4 cases.
1400  double entries = energyResolution[ireg][ie][nResType - 1]->Integral();
1401  for (int ires = 0; ires < nResType; ires++) {
1402  TH1F* hEnergy = (TH1F*)energyResolution[ireg][ie][ires];
1403  double peak = -1.;
1404  double res68 = 99.; // resolution based on 68.3% of the entries
1405  int nIter = 0; // keep track of attempts to fit this histogram
1406  bool fitHist = entries > minEntries;
1407 
1408  //..Possibly iterate fit starting from this point
1409  while (fitHist) {
1410 
1411  //..Fit parameters
1412  double norm = hEnergy->GetMaximum();
1413  double target = fracEnt[nIter] * entries; // fit range contains 68.3% or 50%
1414  std::vector<double> startingParameters;// peak, sigma, eta, fitLow, fitHigh, nbins
1415  startingParameters = eclLeakageFitParameters(hEnergy, target);
1416  func->SetParameters(norm, startingParameters[0], startingParameters[1], startingParameters[2]);
1417  func->SetParLimits(1, startingParameters[3], startingParameters[4]);
1418  func->SetParLimits(2, 0., startingParameters[4] - startingParameters[3]);
1419  func->SetParLimits(3, etaMin, etaMax);
1420 
1421  //..First iteration the fitting range contains >68.3% of the histogram integral.
1422  if (nIter == 0) {
1423 
1424  //..Find the bin numbers used in the fit
1425  int minLo = hEnergy->GetXaxis()->FindBin(startingParameters[3] + 0.0001);
1426  int minHi = hEnergy->GetXaxis()->FindBin(startingParameters[4] - 0.0001);
1427 
1428  //..Subtract a partial bin to get to exactly 68.3%
1429  double dx = hEnergy->GetBinLowEdge(minLo + 1) - hEnergy->GetBinLowEdge(minLo);
1430  double overage = hEnergy->Integral(minLo, minHi) - target;
1431  double subLo = overage / hEnergy->GetBinContent(minLo);
1432  double subHi = overage / hEnergy->GetBinContent(minHi);
1433 
1434  //..Resolution is half the range that contains 68.3% of the events
1435  res68 = 0.5 * dx * (1 + minHi - minLo - max(subLo, subHi));
1436  }
1437 
1438  //..Fit
1439  name = hEnergy->GetName();
1440  B2DEBUG(40, "Fitting " << name.Data());
1441  hEnergy->Fit(func, "LIEQ", "", startingParameters[3], startingParameters[4]);
1442  peak = func->GetParameter(1);
1443  double effSigma = func->GetParameter(2);
1444  double eta = func->GetParameter(3);
1445  double prob = func->GetProb();
1446 
1447  //..Check fit quality 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb,
1448  // 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit
1449  double dEta = min((etaMax - eta), (eta - etaMin));
1450  int fitStatus = eclLeakageFitQuality(startingParameters[3], startingParameters[4], peak, effSigma, dEta, prob);
1451 
1452  //..If the fit probability is low, refit using a smaller range (fracEnt)
1453  if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
1454  nIter++;
1455  } else {
1456  fitHist = false;
1457  }
1458 
1459  } // end of while
1460 
1461  //..Record peak and resolution
1462  peakEnergy[ireg][ie][ires] = peak;
1463  energyRes[ireg][ie][ires] = res68;
1464 
1465  //..And write to disk
1466  histfile->cd();
1467  hEnergy->Write();
1468  }
1469  }
1470  }
1471 
1472  //-----------------------------------------------------------------------------------
1473  //..Summarize resolution
1474  int nresBins = nEnergies * nLeakReg * (nResType + 1); // +1 to add an empty bin after each set
1475  TH1F* peakSummary = new TH1F("peakSummary", "Peak E/Etrue for each method, region, energy;Energy energy point", nresBins, 0,
1476  nEnergies);
1477  TH1F* resolutionSummary = new TH1F("resolutionSummary", "Resolution/peak for each method, region, energy;Energy energy point",
1478  nresBins, 0, nEnergies);
1479 
1480  B2INFO("Resolution divided by peak for each energy bin and region " << nResType << " ways");
1481  for (int ires = 0; ires < nResType; ires++) {B2INFO(" " << resName[ires]);}
1482  ix = 0;
1483  for (int ie = 0; ie < nEnergies; ie++) {
1484  B2INFO(" energy point " << ie);
1485  for (int ireg = 0; ireg < nLeakReg; ireg++) {
1486  B2INFO(" region " << ireg);
1487  for (int ires = 0; ires < nResType; ires++) {
1488 
1489  //..Store in summary histograms
1490  ix++;
1491  peakSummary->SetBinContent(ix, peakEnergy[ireg][ie][ires]);
1492  peakSummary->SetBinError(ix, 0.);
1493  resolutionSummary->SetBinContent(ix, energyRes[ireg][ie][ires] / peakEnergy[ireg][ie][ires]);
1494  resolutionSummary->SetBinError(ix, 0.);
1495 
1496  //..Print out as well
1497  B2INFO(" " << ires << " " << energyRes[ireg][ie][ires] / peakEnergy[ireg][ie][ires]);
1498  }
1499  ix++;
1500  }
1501  }
1502 
1503  //..Write the summary histograms to disk
1504  histfile->cd();
1505  peakSummary->Write();
1506  resolutionSummary->Write();
1507 
1508 
1509  //====================================================================================
1510  //====================================================================================
1511  //..Step 9. Finish up.
1512 
1513  //-----------------------------------------------------------------------------------
1514  //..Close histogram file
1515  histfile->Close();
1516 
1517  //------------------------------------------------------------------------
1518  //..Create and store payload
1519  ECLLeakageCorrections* leakagePayload = new ECLLeakageCorrections();
1520  leakagePayload->setlogEnergiesFwd(forwardVector);
1521  leakagePayload->setlogEnergiesBrl(barrelVector);
1522  leakagePayload->setlogEnergiesBwd(backwardVector);
1523  leakagePayload->setThetaCorrections(thetaCorrection);
1524  leakagePayload->setPhiCorrections(phiCorrection);
1525  leakagePayload->setnCrystalCorrections(nCrystalCorrection);
1526  saveCalibration(leakagePayload, "ECLLeakageCorrections");
1527 
1528  B2INFO("eclLeakageAlgorithm: successfully stored payload ECLLeakageCorrections");
1529  return c_OK;
1530 }
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:32
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:52
void setInitializeActive(bool active)
Setter for m_initializeActive.
Definition: DataStore.cc:92
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 (without leakage correction) divided by generated
float t_locationError
reconstructed minus generated position (cm)
int t_energyBin
generated energy point
double m_noNCrysThreshold
no nCrys fits below this value
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
measured energy with leakage correction 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:95
bool construct(Args &&... params)
Construct an object of type T in this StoreObjPtr, using the provided constructor arguments.
Definition: StoreObjPtr.h:118
static DBStore & Instance()
Instance of a singleton DBStore.
Definition: DBStore.cc:26
void update()
Updates all objects that are outside their interval of validity.
Definition: DBStore.cc:77
Abstract base class for different kinds of events.