Belle II Software  release-05-01-25
eclGammaGammaEAlgorithm.cc
1 
2 #include <ecl/calibration/eclGammaGammaEAlgorithm.h>
3 #include <ecl/dbobjects/ECLCrystalCalib.h>
4 
5 #include "TH2F.h"
6 #include "TFile.h"
7 #include "TMath.h"
8 #include "TF1.h"
9 #include "TROOT.h"
10 
11 using namespace std;
12 using namespace Belle2;
13 using namespace ECL;
14 
17 double eclGammaGammaNovoConst(double* x, double* par)
18 {
19  double qc = 0;
20 
21  double peak = par[1];
22  double width = par[2];
23  double sln4 = sqrt(log(4));
24  double y = par[3] * sln4;
25  double tail = -log(y + sqrt(1 + y * y)) / sln4;
26 
27  if (TMath::Abs(tail) < 1.e-7) {
28  qc = 0.5 * TMath::Power(((x[0] - peak) / width), 2);
29  } else {
30  double qa = tail * sqrt(log(4.));
31  double qb = sinh(qa) / qa;
32  double qx = (x[0] - peak) / width * qb;
33  double qy = 1. + tail * qx;
34 
35  if (qy > 1.E-7)
36  qc = 0.5 * (TMath::Power((log(qy) / tail), 2) + tail * tail);
37  else
38  qc = 15.0;
39  }
40  return par[0] * exp(-qc) + par[4];
41 }
42 
43 
45 eclGammaGammaEAlgorithm::eclGammaGammaEAlgorithm(): CalibrationAlgorithm("eclGammaGammaECollector")
46 {
48  "Perform energy calibration of ecl crystals by fitting a Novosibirsk function to energy deposited by photons in e+e- --> gamma gamma"
49  );
50 }
51 
53 {
56  const double limitTol = 0.0005; /*< tolerance for checking if a parameter is at the limit */
57  const double minFitLimit = 1e-25; /*< cut off for labeling a fit as poor */
58  const double minFitProbIter = 1e-8; /*< cut off for labeling a fit as poor if it also has many iterations */
59  const double constRatio = 0.5; /*< Novosibirsk normalization must be greater than constRatio x constant term */
60  const double peakMin(0.4), peakMax(2.2); /*< range for peak of normalized energy distribution */
61  const double peakTol = limitTol * (peakMax - peakMin); /*< fit is at limit if it is within peakTol of min or max */
62  const double effSigMin(0.02), effSigMax(0.4); /*< range for effective sigma of normalized energy distribution */
63  const double effSigTol = limitTol * (effSigMax - effSigMin); /*< fit is at limit if it is within effSigTol of min or max */
64  const double etaNom(0.41); /*< Nominal tail parameter */
65  const double etaMin(0.), etaMax(5.0); /*< Novosibirsk tail parameter range */
66  const double etaTol = limitTol * (etaMax - etaMin); /*< fit is at limit if it is within etaTol of min or max */
67  const double constTol = 0.001; /*< if constant is less than constTol, it will be fixed to 0 */
68 
70  gROOT->SetBatch();
71 
74  B2INFO("eclGammaGammaAlgorithm parameters:");
75  B2INFO("outputName = " << m_outputName);
76  B2INFO("cellIDLo = " << m_cellIDLo);
77  B2INFO("cellIDHi = " << m_cellIDHi);
78  B2INFO("minEntries = " << m_minEntries);
79  B2INFO("maxIterations = " << m_maxIterations);
80  B2INFO("tRatioMin = " << m_tRatioMin);
81  B2INFO("tRatioMax = " << m_tRatioMax);
82  B2INFO("upperEdgeThresh = " << m_upperEdgeThresh);
83  B2INFO("performFits = " << m_performFits);
84  B2INFO("findExpValues = " << m_findExpValues);
85  B2INFO("storeConst = " << m_storeConst);
86 
89  TH1F* dummy;
90  dummy = (TH1F*)gROOT->FindObject("IntegralVsCrysID");
91  if (dummy) {delete dummy;}
92  dummy = (TH1F*)gROOT->FindObject("AverageExpECrys");
93  if (dummy) {delete dummy;}
94  dummy = (TH1F*)gROOT->FindObject("AverageElecCalib");
95  if (dummy) {delete dummy;}
96  dummy = (TH1F*)gROOT->FindObject("AverageInitCalib");
97  if (dummy) {delete dummy;}
98 
101  auto EnVsCrysID = getObjectPtr<TH2F>("EnVsCrysID");
102  auto ExpEvsCrys = getObjectPtr<TH1F>("ExpEvsCrys");
103  auto ElecCalibvsCrys = getObjectPtr<TH1F>("ElecCalibvsCrys");
104  auto InitialCalibvsCrys = getObjectPtr<TH1F>("InitialCalibvsCrys");
105  auto CalibEntriesvsCrys = getObjectPtr<TH1F>("CalibEntriesvsCrys");
106 
110  TH1F* IntegralVsCrysID = new TH1F("IntegralVsCrysID", "Integral of EnVsCrysID for each crystal;crystal ID;Entries", 8736, 0, 8736);
111  TH1F* AverageExpECrys = new TH1F("AverageExpECrys", "Average expected E per crys from collector;Crystal ID;Energy (GeV)", 8736, 0,
112  8736);
113  TH1F* AverageElecCalib = new TH1F("AverageElecCalib", "Average electronics calib const vs crystal;Crystal ID;Calibration constant",
114  8736, 0, 8736);
115  TH1F* AverageInitCalib = new TH1F("AverageInitCalib", "Average initial calib const vs crystal;Crystal ID;Calibration constant",
116  8736, 0, 8736);
117 
118  for (int crysID = 0; crysID < 8736; crysID++) {
119  TH1D* hEnergy = EnVsCrysID->ProjectionY("hEnergy", crysID + 1, crysID + 1);
120  int Integral = hEnergy->Integral();
121  IntegralVsCrysID->SetBinContent(crysID + 1, Integral);
122 
123  double TotEntries = CalibEntriesvsCrys->GetBinContent(crysID + 1);
124 
125  double expectedE = 0.;
126  if (TotEntries > 0.) {expectedE = ExpEvsCrys->GetBinContent(crysID + 1) / TotEntries;}
127  AverageExpECrys->SetBinContent(crysID + 1, expectedE);
128 
129  double calibconst = 0.;
130  if (TotEntries > 0.) {calibconst = ElecCalibvsCrys->GetBinContent(crysID + 1) / TotEntries;}
131  AverageElecCalib->SetBinContent(crysID + 1, calibconst);
132 
133  calibconst = 0.;
134  if (TotEntries > 0.) {calibconst = InitialCalibvsCrys->GetBinContent(crysID + 1) / TotEntries;}
135  AverageInitCalib->SetBinContent(crysID + 1, calibconst);
136  }
137 
140  TString fName = m_outputName;
141  TFile* histfile = new TFile(fName, "recreate");
142  EnVsCrysID->Write();
143  IntegralVsCrysID->Write();
144  AverageExpECrys->Write();
145  AverageElecCalib->Write();
146  AverageInitCalib->Write();
147 
150  if (!m_performFits) {
151  B2INFO("eclGammaGammaEAlgorithm has not been asked to perform fits; copying input histograms and quitting");
152  histfile->Close();
153  return c_NotEnoughData;
154  }
155 
158  bool sufficientData = true;
159  for (int crysID = m_cellIDLo - 1; crysID < m_cellIDHi; crysID++) {
160  if (IntegralVsCrysID->GetBinContent(crysID + 1) < m_minEntries) {
161  if (m_storeConst == 1) {B2INFO("eclGammaGammaEAlgorithm: crystal " << crysID << " has insufficient statistics: " << IntegralVsCrysID->GetBinContent(crysID + 1) << ". Requirement is " << m_minEntries);}
162  sufficientData = false;
163  break;
164  }
165  }
166 
168  if (!sufficientData && m_storeConst == 1) {
169  histfile->Close();
170  return c_NotEnoughData;
171  }
172 
177  TH1F* CalibVsCrysID = new TH1F("CalibVsCrysID", "Calibration constant vs crystal ID;crystal ID;counts per GeV", 8736, 0, 8736);
178  TH1F* ExpEnergyperCrys = new TH1F("ExpEnergyperCrys", "Expected energy per crystal;Crystal ID;Peak energy (GeV)", 8736, 0, 8736);
179 
181  TH1F* PeakVsCrysID = new TH1F("PeakVsCrysID", "Peak of Novo fit vs crystal ID;crystal ID;Peak normalized energy", 8736, 0, 8736);
182  TH1F* EdgeVsCrysID = new TH1F("EdgeVsCrysID", "Upper edge of Novo fit vs crystal ID;crystal ID;Maximum normalized energy", 8736, 0,
183  8736);
184  TH1F* effSigVsCrysID = new TH1F("effSigVsCrysID", "effSigma vs crystal ID;crystal ID;sigma)", 8736, 0, 8736);
185  TH1F* etaVsCrysID = new TH1F("etaVsCrysID", "eta vs crystal ID;crystal ID;Novo eta parameter", 8736, 0, 8736);
186  TH1F* constVsCrysID = new TH1F("constVsCrysID", "fit constant vs crystal ID;crystal ID;fit constant", 8736, 0, 8736);
187  TH1F* normVsCrysID = new TH1F("normVsCrysID", "Novosibirsk normalization vs crystal ID;crystal ID;normalization", 8736, 0, 8736);
188  TH1F* fitLimitVsCrysID = new TH1F("fitLimitVsCrysID", "fit range lower limit vs crystal ID;crystal ID;upper fit limit", 8736, 0,
189  8736);
190  TH1F* StatusVsCrysID = new TH1F("StatusVsCrysID", "Fit status vs crystal ID;crystal ID;Fit status", 8736, 0, 8736);
191  TH1F* FitProbVsCrysID = new TH1F("FitProbVsCrysID", "Fit probability vs crystal id;crystal ID;Fit probability", 8736, 0, 8736);
192 
194  TH1F* hStatus = new TH1F("hStatus", "Fit status", 25, -5, 20);
195  TH1F* hPeak = new TH1F("hPeak", "Peaks of normalized energy distributions, successful fits;Peak of Novosibirsk fit", 200, 0.8, 1.2);
196  TH1F* fracPeakUnc = new TH1F("fracPeakUnc", "Fractional uncertainty on peak uncertainty, successful fits;Uncertainty on peak", 100,
197  0, 0.1);
198  TH1F* nIterations = new TH1F("nIterations", "Number of times histogram was fit;Number of iterations", 20, -0.5, 19.5);
199 
200 
203  bool allFitsOK = true;
204  for (int crysID = m_cellIDLo - 1; crysID < m_cellIDHi; crysID++) {
205 
207  TString name = "Enormalized";
208  name += crysID;
209  TH1D* hEnergy = EnVsCrysID->ProjectionY(name, crysID + 1, crysID + 1);
210 
212  double histMin = hEnergy->GetXaxis()->GetXmin();
213  double histMax = hEnergy->GetXaxis()->GetXmax();
214  TF1* func = new TF1("eclGammaGammaNovoConst", eclGammaGammaNovoConst, histMin, histMax, 5);
215  func->SetParNames("normalization", "peak", "effSigma", "eta", "const");
216  func->SetParLimits(1, peakMin, peakMax);
217  func->SetParLimits(2, effSigMin, effSigMax);
218  func->SetParLimits(3, etaMin, etaMax);
219  func->SetParLimits(4, 0., hEnergy->GetMaximum());
220 
222  hEnergy->GetXaxis()->SetRangeUser(peakMin, peakMax);
223  int maxBin = hEnergy->GetMaximumBin();
224  double peakE = hEnergy->GetBinLowEdge(maxBin);
225  double peakEUnc = 0.;
226  double normalization = hEnergy->GetMaximum();
227  double normUnc = 0.;
228  double effSigma = hEnergy->GetRMS();
229  double sigmaUnc = 0.;
230  hEnergy->GetXaxis()->SetRangeUser(histMin, histMax);
231 
233  double fitlow = peakE - effSigma;
234  double fithigh = histMax;
235 
237  double eta = etaNom;
238  double etaUnc = 0.;
239  double constant = 0.01 * normalization;
240  double constUnc = 0.;
241 
243  double dIter = 0.1 * (histMax - histMin) / hEnergy->GetNbinsX();
244  double fitProb(0.);
245  double fitProbDefault(0.);
246  double lowold(0.), lowoldold(0.);
247  bool fixConst = false;
248  int nIter = 0;
249  bool fitHist = IntegralVsCrysID->GetBinContent(crysID + 1) >= m_minEntries; /* fit only if enough events */
250 
253  while (fitHist) {
254 
256  func->SetParameters(normalization, peakE, effSigma, eta, constant);
257  if (fixConst) { func->FixParameter(4, 0); }
258 
260  hEnergy->Fit(func, "LIQ", "", fitlow, fithigh);
261  nIter++;
262  fitHist = false;
263  normalization = func->GetParameter(0);
264  normUnc = func->GetParError(0);
265  peakE = func->GetParameter(1);
266  peakEUnc = func->GetParError(1);
267  effSigma = func->GetParameter(2);
268  sigmaUnc = func->GetParError(2);
269  eta = func->GetParameter(3);
270  etaUnc = func->GetParError(3);
271  constant = func->GetParameter(4);
272  constUnc = func->GetParError(4);
273  fitProbDefault = func->GetProb();
274 
276  double peak = func->Eval(peakE) - constant;
277  double tRatio = (func->Eval(fitlow) - constant) / peak;
278  if (tRatio < m_tRatioMin || tRatio > m_tRatioMax) {
279  double targetY = constant + 0.5 * (m_tRatioMin + m_tRatioMax) * peak;
280  lowoldold = lowold;
281  lowold = fitlow;
282  fitlow = func->GetX(targetY, histMin, peakE);
283  fitHist = true;
284 
286  if (abs(fitlow - lowoldold) < dIter) {fitlow = 0.5 * (lowold + lowoldold); }
287 
289  if (nIter > m_maxIterations - 3) {fitlow = 0.33333 * (fitlow + lowold + lowoldold); }
290  }
291 
293  if (constant < constTol && !fixConst) {
294  constant = 0;
295  fixConst = true;
296  fitHist = true;
297  }
298 
300  if (nIter == m_maxIterations) {fitHist = false;}
301  B2DEBUG(10, crysID << " " << nIter << " " << peakE << " " << constant << " " << tRatio << " " << fitlow);
302  }
303 
306  fitProb = 0.;
307  if (nIter > 0) {
308  int lowbin = hEnergy->GetXaxis()->FindBin(fitlow);
309  int highbin = hEnergy->GetXaxis()->FindBin(fithigh);
310  int npar = 5;
311  if (fixConst) {npar = 4;}
312  int ndeg = -npar;
313  double chisq = 0.;
314  double binwidth = hEnergy->GetBinWidth(1);
315  for (int ib = lowbin; ib <= highbin; ib++) {
316  double xlow = hEnergy->GetBinLowEdge(ib);
317  double yexp = func->Integral(xlow, xlow + binwidth) / binwidth;
318 
320  if (yexp > constTol) {
321  double yobs = hEnergy->GetBinContent(ib);
322  double dnom = yexp;
323  if (yexp < 0.9999 && yobs > yexp) {dnom = yobs;}
324  double dchi2 = (yexp - yobs) * (yexp - yobs) / dnom;
325  chisq += dchi2;
326  ndeg++;
327  }
328  }
329  fitProb = TMath::Prob(chisq, ndeg);
330  }
331 
334  int iStatus = fitOK; // success
335  if (nIter == m_maxIterations) {iStatus = iterations;} // too many iterations
336 
338  if (normalization < constRatio * constant) {iStatus = noPeak;}
339 
341  if (fitProb <= minFitLimit || (fitProb < minFitProbIter && iStatus == iterations)) {iStatus = poorFit;}
342 
344  if ((peakE < peakMin + peakTol) || (peakE > peakMax - peakTol)) {iStatus = atLimit;}
345  if ((effSigma < effSigMin + effSigTol) || (effSigma > effSigMax - effSigTol)) {iStatus = atLimit;}
346  if ((eta < etaMin + etaTol) || (eta > etaMax - etaTol)) {iStatus = atLimit;}
347 
348  //** No fit
349  if (nIter == 0) {iStatus = notFit;} // not fit
350 
353  double upperEdge = peakE;
354  double edgeUnc = peakEUnc;
355 
356  if (iStatus >= iterations) {
357 
359  double targetY = constant + m_upperEdgeThresh * (func->Eval(peakE) - constant);
360 
362  int iLow = hEnergy->GetXaxis()->FindBin(peakE) + 1;
363  int iHigh = hEnergy->GetNbinsX();
364  int iLast = iLow;
365  for (int ibin = iLow; ibin < iHigh; ibin++) {
366  double xc = hEnergy->GetBinCenter(ibin);
367  if (func->Eval(xc) > targetY) {iLast = ibin;}
368  }
369  double xLow = hEnergy->GetBinCenter(iLast);
370  double xHigh = hEnergy->GetBinCenter(iLast + 1);
371 
373  func->SetNpx(1000);
374  upperEdge = func->GetX(targetY, xLow, xHigh);
375 
376 
377  } else if (iStatus > notFit) {
378 
380  int iLast = -1;
381  int thisBin = hEnergy->GetBinContent(1);
382  for (int ibin = 2; ibin < hEnergy->GetNbinsX(); ibin++) {
383  int prevBin = thisBin;
384  thisBin = hEnergy->GetBinContent(ibin);
385  if (thisBin > 0 && thisBin + prevBin >= 2) {iLast = ibin;}
386  }
387  if (iLast > 0) {upperEdge = hEnergy->GetBinCenter(iLast);}
388  }
389 
392  int histbin = crysID + 1;
393  PeakVsCrysID->SetBinContent(histbin, peakE);
394  PeakVsCrysID->SetBinError(histbin, peakEUnc);
395  EdgeVsCrysID->SetBinContent(histbin, upperEdge);
396  EdgeVsCrysID->SetBinError(histbin, edgeUnc);
397  effSigVsCrysID->SetBinContent(histbin, effSigma);
398  effSigVsCrysID->SetBinError(histbin, sigmaUnc);
399  etaVsCrysID->SetBinContent(histbin, eta);
400  etaVsCrysID->SetBinError(histbin, etaUnc);
401  constVsCrysID->SetBinContent(histbin, constant);
402  constVsCrysID->SetBinError(histbin, constUnc);
403  normVsCrysID->SetBinContent(histbin, normalization);
404  normVsCrysID->SetBinError(histbin, normUnc);
405  fitLimitVsCrysID->SetBinContent(histbin, fitlow);
406  fitLimitVsCrysID->SetBinError(histbin, 0);
407  StatusVsCrysID->SetBinContent(histbin, iStatus);
408  StatusVsCrysID->SetBinError(histbin, 0);
409  FitProbVsCrysID->SetBinContent(histbin, fitProb);
410  FitProbVsCrysID->SetBinError(histbin, 0);
411 
413  hStatus->Fill(iStatus);
414  nIterations->Fill(nIter);
415  if (iStatus >= iterations) {
416  hPeak->Fill(peakE);
417  fracPeakUnc->Fill(peakEUnc / peakE);
418  }
419 
421  B2INFO("cellID " << crysID + 1 << " status = " << iStatus << " fit probability = " << fitProb << " default prob = " <<
422  fitProbDefault);
423  histfile->cd();
424  hEnergy->Write();
425 
426  } /* end of loop over crystals */
427 
430  for (int crysID = 0; crysID < 8736; crysID++) {
431  int histbin = crysID + 1;
432  double fitstatus = StatusVsCrysID->GetBinContent(histbin);
433  double upperEdge = EdgeVsCrysID->GetBinContent(histbin);
434  double fracEdgeUnc = EdgeVsCrysID->GetBinError(histbin) / upperEdge;
435 
437  if (fitstatus < 0) {
438  upperEdge = -1.;
439  fracEdgeUnc = 0.;
440  if (histbin >= m_cellIDLo && histbin <= m_cellIDHi) {
441  B2INFO("eclGammaGammaEAlgorithm: cellID " << histbin << " is not a successful fit. Status = " << fitstatus);
442  allFitsOK = false;
443  }
444  }
445 
447  if (m_findExpValues) {
448  double inputExpE = abs(AverageExpECrys->GetBinContent(histbin));
449  ExpEnergyperCrys->SetBinContent(histbin, inputExpE * upperEdge);
450  ExpEnergyperCrys->SetBinError(histbin, fracEdgeUnc * inputExpE * upperEdge);
451  } else {
452 
454  double inputCalib = abs(AverageInitCalib->GetBinContent(histbin));
455  CalibVsCrysID->SetBinContent(histbin, inputCalib / upperEdge);
456  CalibVsCrysID->SetBinError(histbin, fracEdgeUnc * inputCalib / upperEdge);
457  }
458  }
459 
462  bool DBsuccess = false;
463  if (m_storeConst == 0 || (m_storeConst == 1 && allFitsOK)) {
464  DBsuccess = true;
465  if (m_findExpValues) {
466 
468  std::vector<float> tempE;
469  std::vector<float> tempUnc;
470  for (int crysID = 0; crysID < 8736; crysID++) {
471  tempE.push_back(ExpEnergyperCrys->GetBinContent(crysID + 1));
472  tempUnc.push_back(ExpEnergyperCrys->GetBinError(crysID + 1));
473  }
474  ECLCrystalCalib* ExpectedE = new ECLCrystalCalib();
475  ExpectedE->setCalibVector(tempE, tempUnc);
476  saveCalibration(ExpectedE, "ECLExpGammaGammaE");
477  B2INFO("eclCosmicEAlgorithm: successfully stored expected energies ECLExpGammaGammaE");
478 
479  } else {
480 
482  std::vector<float> tempCalib;
483  std::vector<float> tempCalibUnc;
484  for (int crysID = 0; crysID < 8736; crysID++) {
485  tempCalib.push_back(CalibVsCrysID->GetBinContent(crysID + 1));
486  tempCalibUnc.push_back(CalibVsCrysID->GetBinError(crysID + 1));
487  }
488  ECLCrystalCalib* GammaGammaECalib = new ECLCrystalCalib();
489  GammaGammaECalib->setCalibVector(tempCalib, tempCalibUnc);
490  saveCalibration(GammaGammaECalib, "ECLCrystalEnergyGammaGamma");
491  B2INFO("eclGammaGammaEAlgorithm: successfully stored ECLCrystalEnergyGammaGamma calibration constants");
492  }
493  }
494 
498  PeakVsCrysID->Write();
499  EdgeVsCrysID->Write();
500  effSigVsCrysID->Write();
501  etaVsCrysID->Write();
502  constVsCrysID->Write();
503  normVsCrysID->Write();
504  fitLimitVsCrysID->Write();
505  StatusVsCrysID->Write();
506  FitProbVsCrysID->Write();
507  hPeak->Write();
508  fracPeakUnc->Write();
509  nIterations->Write();
510  hStatus->Write();
511 
513  if (m_findExpValues) {
514  ExpEnergyperCrys->Write();
515  } else {
516  CalibVsCrysID->Write();
517  }
518  histfile->Close();
519 
522  dummy = (TH1F*)gROOT->FindObject("PeakVsCrysID"); delete dummy;
523  dummy = (TH1F*)gROOT->FindObject("EdgeVsCrysID"); delete dummy;
524  dummy = (TH1F*)gROOT->FindObject("effSigVsCrysID"); delete dummy;
525  dummy = (TH1F*)gROOT->FindObject("etaVsCrysID"); delete dummy;
526  dummy = (TH1F*)gROOT->FindObject("constVsCrysID"); delete dummy;
527  dummy = (TH1F*)gROOT->FindObject("normVsCrysID"); delete dummy;
528  dummy = (TH1F*)gROOT->FindObject("fitLimitVsCrysID"); delete dummy;
529  dummy = (TH1F*)gROOT->FindObject("StatusVsCrysID"); delete dummy;
530  dummy = (TH1F*)gROOT->FindObject("FitProbVsCrysID"); delete dummy;
531  dummy = (TH1F*)gROOT->FindObject("fracPeakUnc"); delete dummy;
532  dummy = (TH1F*)gROOT->FindObject("nIterations"); delete dummy;
533  dummy = (TH1F*)gROOT->FindObject("hStatus"); delete dummy;
534  dummy = (TH1F*)gROOT->FindObject("ExpEnergyperCrys"); delete dummy;
535  dummy = (TH1F*)gROOT->FindObject("CalibVsCrysID"); delete dummy;
536 
537 
540  if (m_storeConst == -1) {
541  B2INFO("eclGammaGammaEAlgorithm performed fits but was not asked to store contants");
542  return c_Failure;
543  } else if (!DBsuccess) {
544  if (m_findExpValues) { B2INFO("eclGammaGammaEAlgorithm: failed to store expected values"); }
545  else { B2INFO("eclGammaGammaEAlgorithm: failed to store calibration constants"); }
546  return c_Failure;
547  }
548  return c_OK;
549 }
Belle2::ECL::eclGammaGammaEAlgorithm::notFit
int notFit
no fit performed; no constants found for this crystal
Definition: eclGammaGammaEAlgorithm.h:145
Belle2::ECL::eclGammaGammaEAlgorithm::m_maxIterations
int m_maxIterations
no more than maxIteration iterations
Definition: eclGammaGammaEAlgorithm.h:127
Belle2::ECL::eclGammaGammaEAlgorithm::m_upperEdgeThresh
double m_upperEdgeThresh
Upper edge is where the fit = upperEdgeThresh * peak value.
Definition: eclGammaGammaEAlgorithm.h:130
Belle2::ECL::eclGammaGammaEAlgorithm::m_cellIDHi
int m_cellIDHi
Last cellID to be fit.
Definition: eclGammaGammaEAlgorithm.h:125
Belle2::ECL::eclGammaGammaEAlgorithm::noPeak
int noPeak
Novosibirsk component of fit is negligible; upper edge is found from histogram, not fit.
Definition: eclGammaGammaEAlgorithm.h:144
Belle2::ECL::eclGammaGammaEAlgorithm::atLimit
int atLimit
a parameter is at the limit; upper edge is found from histogram, not fit
Definition: eclGammaGammaEAlgorithm.h:142
Belle2::ECL::eclGammaGammaEAlgorithm::fitOK
int fitOK
Characterize fit status.
Definition: eclGammaGammaEAlgorithm.h:140
Belle2::ECL::eclGammaGammaEAlgorithm::iterations
int iterations
fit reached max number of iterations, but is useable
Definition: eclGammaGammaEAlgorithm.h:141
Belle2::CalibrationAlgorithm::saveCalibration
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
Definition: CalibrationAlgorithm.cc:290
Belle2::ECL::eclGammaGammaEAlgorithm::m_performFits
bool m_performFits
if false, input histograms are copied to output, but no fits are done
Definition: eclGammaGammaEAlgorithm.h:131
Belle2::ECL::eclGammaGammaEAlgorithm::calibrate
virtual EResult calibrate() override
..Run algorithm on events
Definition: eclGammaGammaEAlgorithm.cc:52
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::ECLCrystalCalib
General DB object to store one calibration number per ECL crystal.
Definition: ECLCrystalCalib.h:34
Belle2::CalibrationAlgorithm::setDescription
void setDescription(const std::string &description)
Set algorithm description (in constructor)
Definition: CalibrationAlgorithm.h:331
Belle2::ECL::eclGammaGammaEAlgorithm::m_tRatioMax
double m_tRatioMax
Fit range is adjusted so that fit at upper endpoint is between tRatioMin and tRatioMax of peak.
Definition: eclGammaGammaEAlgorithm.h:129
Belle2::ECL::eclGammaGammaEAlgorithm::m_minEntries
int m_minEntries
Minimum entries to fit a crystal.
Definition: eclGammaGammaEAlgorithm.h:126
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::ECL::eclGammaGammaEAlgorithm::m_findExpValues
bool m_findExpValues
if true, fits are used to find expected energy deposit for each crystal instead of the calibration co...
Definition: eclGammaGammaEAlgorithm.h:132
Belle2::CalibrationAlgorithm::c_Failure
@ c_Failure
Failed =3 in Python.
Definition: CalibrationAlgorithm.h:54
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::ECL::eclGammaGammaEAlgorithm::m_outputName
std::string m_outputName
..Parameters to control Novosibirsk fit to energy deposited in each crystal by mu+mu- events
Definition: eclGammaGammaEAlgorithm.h:123
Belle2::CalibrationAlgorithm::c_NotEnoughData
@ c_NotEnoughData
Needs more data =2 in Python.
Definition: CalibrationAlgorithm.h:53
Belle2::ECLCrystalCalib::setCalibVector
void setCalibVector(const std::vector< float > &CalibConst, const std::vector< float > &CalibConstUnc)
Set vector of constants with uncertainties.
Definition: ECLCrystalCalib.h:48
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47
Belle2::ECL::eclGammaGammaEAlgorithm::m_tRatioMin
double m_tRatioMin
Fit range is adjusted so that fit at upper endpoint is between tRatioMin and tRatioMax of peak.
Definition: eclGammaGammaEAlgorithm.h:128
Belle2::ECL::eclGammaGammaEAlgorithm::m_cellIDLo
int m_cellIDLo
First cellID to be fit.
Definition: eclGammaGammaEAlgorithm.h:124
Belle2::ECL::eclGammaGammaEAlgorithm::poorFit
int poorFit
low chi square; upper edge is found from histogram, not fit
Definition: eclGammaGammaEAlgorithm.h:143
Belle2::ECL::eclGammaGammaEAlgorithm::m_storeConst
int m_storeConst
controls which values are written to the database.
Definition: eclGammaGammaEAlgorithm.h:134