Belle II Software  release-05-01-25
eclMuMuEAlgorithm.cc
1 
2 #include <ecl/calibration/eclMuMuEAlgorithm.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 
16 double eclNovoConst(double* x, double* par)
17 {
18  double qc = 0;
19 
20  double peak = par[1];
21  double width = par[2];
22  double sln4 = sqrt(log(4));
23  double y = par[3] * sln4;
24  double tail = -log(y + sqrt(1 + y * y)) / sln4;
25 
26  if (TMath::Abs(tail) < 1.e-7) {
27  qc = 0.5 * TMath::Power(((x[0] - peak) / width), 2);
28  } else {
29  double qa = tail * sqrt(log(4.));
30  double qb = sinh(qa) / qa;
31  double qx = (x[0] - peak) / width * qb;
32  double qy = 1. + tail * qx;
33 
34  if (qy > 1.E-7)
35  qc = 0.5 * (TMath::Power((log(qy) / tail), 2) + tail * tail);
36  else
37  qc = 15.0;
38  }
39  return par[0] * exp(-qc) + par[4];
40 }
41 
42 
43 eclMuMuEAlgorithm::eclMuMuEAlgorithm(): CalibrationAlgorithm("eclMuMuECollector"), cellIDLo(1), cellIDHi(8736), minEntries(150),
44  maxIterations(10), tRatioMin(0.2), tRatioMax(0.25), performFits(true), findExpValues(false), storeConst(0)
45 {
47  "Perform energy calibration of ecl crystals by fitting a Novosibirsk function to energy deposited by muons"
48  );
49 }
50 
52 {
55  double limitTol = 0.0005; /*< tolerance for checking if a parameter is at the limit */
56  double minFitLimit = 1e-25; /*< cut off for labeling a fit as poor */
57  double minFitProbIter = 1e-8; /*< cut off for labeling a fit as poor if it also has many iterations */
58  double constRatio = 0.5; /*< Novosibirsk normalization must be greater than constRatio x constant term */
59  double peakMin(0.4), peakMax(2.2); /*< range for peak of normalized energy distribution */
60  double peakTol = limitTol * (peakMax - peakMin); /*< fit is at limit if it is within peakTol of min or max */
61  double effSigMin(0.02), effSigMax(0.2); /*< range for effective sigma of normalized energy distribution */
62  double effSigTol = limitTol * (effSigMax - effSigMin); /*< fit is at limit if it is within effSigTol of min or max */
63  double etaNom(-0.41); /*< Nominal tail parameter; fixed to this value for endcap and when findExValues=false */
64  double etaMin(-0.7), etaMax(-0.2); /*< Novosibirsk tail parameter range */
65  double etaTol = limitTol * (etaMax - etaMin); /*< fit is at limit if it is within etaTol of min or max */
66  double constMin(0.), constMax(10.); /*< constant term in normalized energy distribution */
67  double constTol = limitTol * constMax; /*< if constant is less than constTol, it will be fixed to 0 */
68 
70  gROOT->SetBatch();
71 
74  TH1F* dummy;
75  dummy = (TH1F*)gROOT->FindObject("IntegralVsCrysID");
76  if (dummy) {delete dummy;}
77  dummy = (TH1F*)gROOT->FindObject("AverageExpECrys");
78  if (dummy) {delete dummy;}
79  dummy = (TH1F*)gROOT->FindObject("AverageElecCalib");
80  if (dummy) {delete dummy;}
81  dummy = (TH1F*)gROOT->FindObject("AverageInitCalib");
82  if (dummy) {delete dummy;}
83 
86  auto EnVsCrysID = getObjectPtr<TH2F>("EnVsCrysID");
87  auto ExpEvsCrys = getObjectPtr<TH1F>("ExpEvsCrys");
88  auto ElecCalibvsCrys = getObjectPtr<TH1F>("ElecCalibvsCrys");
89  auto InitialCalibvsCrys = getObjectPtr<TH1F>("InitialCalibvsCrys");
90  auto CalibEntriesvsCrys = getObjectPtr<TH1F>("CalibEntriesvsCrys");
91  auto RawDigitAmpvsCrys = getObjectPtr<TH2F>("RawDigitAmpvsCrys");
92  auto RawDigitTimevsCrys = getObjectPtr<TH2F>("RawDigitTimevsCrys");
93 
97  TH1F* IntegralVsCrysID = new TH1F("IntegralVsCrysID", "Integral of EnVsCrysID for each crystal;crystal ID;Entries", 8736, 0, 8736);
98  TH1F* AverageExpECrys = new TH1F("AverageExpECrys", "Average expected E per crys from collector;Crystal ID;Energy (GeV)", 8736, 0,
99  8736);
100  TH1F* AverageElecCalib = new TH1F("AverageElecCalib", "Average electronics calib const vs crystal;Crystal ID;Calibration constant",
101  8736, 0, 8736);
102  TH1F* AverageInitCalib = new TH1F("AverageInitCalib", "Average initial calib const vs crystal;Crystal ID;Calibration constant",
103  8736, 0, 8736);
104 
105  for (int crysID = 0; crysID < 8736; crysID++) {
106  TH1D* hEnergy = EnVsCrysID->ProjectionY("hEnergy", crysID + 1, crysID + 1);
107  int Integral = hEnergy->Integral();
108  IntegralVsCrysID->SetBinContent(crysID + 1, Integral);
109 
110  double TotEntries = CalibEntriesvsCrys->GetBinContent(crysID + 1);
111 
112  double expectedE = 0.;
113  if (TotEntries > 0.) {expectedE = ExpEvsCrys->GetBinContent(crysID + 1) / TotEntries;}
114  AverageExpECrys->SetBinContent(crysID + 1, expectedE);
115 
116  double calibconst = 0.;
117  if (TotEntries > 0.) {calibconst = ElecCalibvsCrys->GetBinContent(crysID + 1) / TotEntries;}
118  AverageElecCalib->SetBinContent(crysID + 1, calibconst);
119 
120  calibconst = 0.;
121  if (TotEntries > 0.) {calibconst = InitialCalibvsCrys->GetBinContent(crysID + 1) / TotEntries;}
122  AverageInitCalib->SetBinContent(crysID + 1, calibconst);
123  }
124 
127  TFile* histfile = new TFile("eclMuMuEAlgorithm.root", "recreate");
128  EnVsCrysID->Write();
129  IntegralVsCrysID->Write();
130  AverageExpECrys->Write();
131  AverageElecCalib->Write();
132  AverageInitCalib->Write();
133  RawDigitAmpvsCrys->Write();
134  RawDigitTimevsCrys->Write();
135 
138  if (!performFits) {
139  B2RESULT("eclMuMuEAlgorithm has not been asked to perform fits; copying input histograms and quitting");
140  histfile->Close();
141  return c_NotEnoughData;
142  }
143 
146  bool sufficientData = true;
147  for (int crysID = cellIDLo - 1; crysID < cellIDHi; crysID++) {
148  if (IntegralVsCrysID->GetBinContent(crysID + 1) < minEntries) {
149  if (storeConst == 1) {B2RESULT("eclMuMuEAlgorithm: crystal " << crysID << " has insufficient statistics: " << IntegralVsCrysID->GetBinContent(crysID + 1) << ". Requirement is " << minEntries);}
150  sufficientData = false;
151  break;
152  }
153  }
154 
156  if (!sufficientData && storeConst == 1) {
157  histfile->Close();
158  return c_NotEnoughData;
159  }
160 
165  TH1F* CalibVsCrysID = new TH1F("CalibVsCrysID", "Calibration constant vs crystal ID;crystal ID;counts per GeV", 8736, 0, 8736);
166  TH1F* ExpEnergyperCrys = new TH1F("ExpEnergyperCrys", "Expected energy per crystal;Crystal ID;Peak energy (GeV)", 8736, 0, 8736);
167 
169  TH1F* PeakVsCrysID = new TH1F("PeakVsCrysID", "Peak of Novo fit vs crystal ID;crystal ID;Peak normalized energy", 8736, 0,
170  8736);
171  TH1F* effSigVsCrysID = new TH1F("effSigVsCrysID", "effSigma vs crystal ID;crystal ID;sigma)", 8736, 0, 8736);
172  TH1F* etaVsCrysID = new TH1F("etaVsCrysID", "eta vs crystal ID;crystal ID;Novo eta parameter", 8736, 0, 8736);
173  TH1F* constVsCrysID = new TH1F("constVsCrysID", "fit constant vs crystal ID;crystal ID;fit constant", 8736, 0, 8736);
174  TH1F* normVsCrysID = new TH1F("normVsCrysID", "Novosibirsk normalization vs crystal ID;crystal ID;normalization", 8736, 0, 8736);
175  TH1F* fitLimitVsCrysID = new TH1F("fitLimitVsCrysID", "fit range upper limit vs crystal ID;crystal ID;upper fit limit", 8736, 0,
176  8736);
177  TH1F* StatusVsCrysID = new TH1F("StatusVsCrysID", "Fit status vs crystal ID;crystal ID;Fit status", 8736, 0, 8736);
178 
180  TH1F* hStatus = new TH1F("hStatus", "Fit status", 25, -5, 20);
181  TH1F* hPeak = new TH1F("hPeak", "Peaks of normalized energy distributions, successful fits;Peak of Novosibirsk fit", 200, 0.8, 1.2);
182  TH1F* fracPeakUnc = new TH1F("fracPeakUnc", "Fractional uncertainty on peak uncertainty, successful fits;Uncertainty on peak", 100,
183  0, 0.1);
184  TH1F* nIterations = new TH1F("nIterations", "Number of times histogram was fit;Number of iterations", 20, -0.5, 19.5);
185 
186 
189  bool allFitsOK = true;
190  for (int crysID = cellIDLo - 1; crysID < cellIDHi; crysID++) {
191 
193  TString name = "Enormalized";
194  name += crysID;
195  TH1D* hEnergy = EnVsCrysID->ProjectionY(name, crysID + 1, crysID + 1);
196 
198  double histMin = hEnergy->GetXaxis()->GetXmin();
199  double histMax = hEnergy->GetXaxis()->GetXmax();
200  TF1* func = new TF1("eclNovoConst", eclNovoConst, histMin, histMax, 5);
201  func->SetParNames("normalization", "peak", "effSigma", "eta", "const");
202  func->SetParLimits(1, peakMin, peakMax);
203  func->SetParLimits(2, effSigMin, effSigMax);
204  func->SetParLimits(3, etaMin, etaMax);
205  func->SetParLimits(4, constMin, constMax);
206 
208  hEnergy->GetXaxis()->SetRangeUser(peakMin, peakMax);
209  int maxBin = hEnergy->GetMaximumBin();
210  double peakE = hEnergy->GetBinLowEdge(maxBin);
211  double peakEUnc = 0.;
212  double normalization = hEnergy->GetMaximum();
213  double normUnc = 0.;
214  double effSigma = hEnergy->GetRMS();
215  double sigmaUnc = 0.;
216  hEnergy->GetXaxis()->SetRangeUser(histMin, histMax);
217 
219  double fitlow = 0.25;
220  double fithigh = peakE + 2.5 * effSigma;
221 
223  double eta = etaNom;
224  double etaUnc = 0.;
225 
227  int il0 = hEnergy->GetXaxis()->FindBin(fitlow);
228  int il1 = hEnergy->GetXaxis()->FindBin(fitlow + 0.1);
229  double constant = hEnergy->Integral(il0, il1) / (1 + il1 - il0);
230  if (constant < 0.01 * normalization) {constant = 0.01 * normalization;}
231  double constUnc = 0.;
232 
234  double dIter = 0.1 * (histMax - histMin) / hEnergy->GetNbinsX();
235  double fitProb(0.);
236  double highold(0.), higholdold(0.);
237  bool fixConst = false;
238  int nIter = 0;
239  bool fitHist = IntegralVsCrysID->GetBinContent(crysID + 1) >= minEntries; /* fit only if enough events */
240 
243  while (fitHist) {
244 
246  func->SetParameters(normalization, peakE, effSigma, eta, constant);
247  if (fixConst) { func->FixParameter(4, 0); }
248 
250  if (crysID < 1152 || crysID > 7775 || !findExpValues) {func->FixParameter(3, etaNom);}
251 
253  hEnergy->Fit(func, "LIQ", "", fitlow, fithigh);
254  nIter++;
255  fitHist = false;
256  normalization = func->GetParameter(0);
257  normUnc = func->GetParError(0);
258  peakE = func->GetParameter(1);
259  peakEUnc = func->GetParError(1);
260  effSigma = func->GetParameter(2);
261  sigmaUnc = func->GetParError(2);
262  eta = func->GetParameter(3);
263  etaUnc = func->GetParError(3);
264  constant = func->GetParameter(4);
265  constUnc = func->GetParError(4);
266  fitProb = func->GetProb();
267 
269  double peak = func->Eval(peakE) - constant;
270  double tRatio = (func->Eval(fithigh) - constant) / peak;
271  if (tRatio < tRatioMin || tRatio > tRatioMax) {
272  double targetY = constant + 0.5 * (tRatioMin + tRatioMax) * peak;
273  higholdold = highold;
274  highold = fithigh;
275  fithigh = func->GetX(targetY, peakE, histMax);
276  fitHist = true;
277 
279  if (abs(fithigh - higholdold) < dIter) {fithigh = 0.5 * (highold + higholdold); }
280 
282  if (nIter > maxIterations - 3) {fithigh = 0.33333 * (fithigh + highold + higholdold); }
283  }
284 
286  if (constant < constTol && !fixConst) {
287  constant = 0;
288  fixConst = true;
289  fitHist = true;
290  }
291 
293  if (nIter == maxIterations) {fitHist = false;}
294  B2DEBUG(200, crysID << " " << nIter << " " << peakE << " " << constant << " " << tRatio << " " << fithigh);
295  }
296 
299  int iStatus = fitOK; // success
300  if (nIter == maxIterations) {iStatus = iterations;} // too many iterations
301 
303  if (normalization < constRatio * constant) {iStatus = noPeak;}
304 
306  if (fitProb <= minFitLimit || (fitProb < minFitProbIter && iStatus == iterations)) {iStatus = poorFit;}
307 
309  if ((peakE < peakMin + peakTol) || (peakE > peakMax - peakTol)) {iStatus = atLimit;}
310  if ((effSigma < effSigMin + effSigTol) || (effSigma > effSigMax - effSigTol)) {iStatus = atLimit;}
311  if ((eta < etaMin + etaTol) || (eta > etaMax - etaTol)) {iStatus = atLimit;}
312  if (constant > constMax - constTol) {iStatus = atLimit;}
313 
314  //** No fit
315  if (nIter == 0) {iStatus = notFit;} // not fit
316 
318  int histbin = crysID + 1;
319  PeakVsCrysID->SetBinContent(histbin, peakE);
320  PeakVsCrysID->SetBinError(histbin, peakEUnc);
321  effSigVsCrysID->SetBinContent(histbin, effSigma);
322  effSigVsCrysID->SetBinError(histbin, sigmaUnc);
323  etaVsCrysID->SetBinContent(histbin, eta);
324  etaVsCrysID->SetBinError(histbin, etaUnc);
325  constVsCrysID->SetBinContent(histbin, constant);
326  constVsCrysID->SetBinError(histbin, constUnc);
327  normVsCrysID->SetBinContent(histbin, normalization);
328  normVsCrysID->SetBinError(histbin, normUnc);
329  fitLimitVsCrysID->SetBinContent(histbin, fithigh);
330  fitLimitVsCrysID->SetBinError(histbin, 0);
331  StatusVsCrysID->SetBinContent(histbin, iStatus);
332  StatusVsCrysID->SetBinError(histbin, 0);
333 
335  hStatus->Fill(iStatus);
336  nIterations->Fill(nIter);
337  if (iStatus >= iterations) {
338  hPeak->Fill(peakE);
339  fracPeakUnc->Fill(peakEUnc / peakE);
340  }
341 
343  B2INFO("cellID " << crysID + 1 << " status = " << iStatus << " fit probability = " << fitProb);
344  histfile->cd();
345  hEnergy->Write();
346 
347  } /* end of loop over crystals */
348 
351  for (int crysID = 0; crysID < 8736; crysID++) {
352  int histbin = crysID + 1;
353  double fitstatus = StatusVsCrysID->GetBinContent(histbin);
354  double peakE = PeakVsCrysID->GetBinContent(histbin);
355  double fracpeakEUnc = PeakVsCrysID->GetBinError(histbin) / peakE;
356 
358  if (fitstatus < iterations) {
359  peakE = -1.;
360  fracpeakEUnc = 0.;
361  if (histbin >= cellIDLo && histbin <= cellIDHi) {
362  B2RESULT("eclMuMuEAlgorithm: cellID " << histbin << " is not a successful fit. Status = " << fitstatus);
363  allFitsOK = false;
364  }
365  }
366 
368  if (findExpValues) {
369  double inputExpE = abs(AverageExpECrys->GetBinContent(histbin));
370  ExpEnergyperCrys->SetBinContent(histbin, inputExpE * peakE);
371  ExpEnergyperCrys->SetBinError(histbin, fracpeakEUnc * inputExpE * peakE);
372  } else {
373 
375  double inputCalib = abs(AverageInitCalib->GetBinContent(histbin));
376  CalibVsCrysID->SetBinContent(histbin, inputCalib / peakE);
377  CalibVsCrysID->SetBinError(histbin, fracpeakEUnc * inputCalib / peakE);
378  }
379  }
380 
383  bool DBsuccess = false;
384  if (storeConst == 0 || (storeConst == 1 && allFitsOK)) {
385  DBsuccess = true;
386  if (findExpValues) {
387 
389  std::vector<float> tempE;
390  std::vector<float> tempUnc;
391  for (int crysID = 0; crysID < 8736; crysID++) {
392  tempE.push_back(ExpEnergyperCrys->GetBinContent(crysID + 1));
393  tempUnc.push_back(ExpEnergyperCrys->GetBinError(crysID + 1));
394  }
395  ECLCrystalCalib* ExpectedE = new ECLCrystalCalib();
396  ExpectedE->setCalibVector(tempE, tempUnc);
397  saveCalibration(ExpectedE, "ECLExpMuMuE");
398  B2RESULT("eclCosmicEAlgorithm: successfully stored expected energies ECLExpMuMuE");
399 
400  } else {
401 
403  std::vector<float> tempCalib;
404  std::vector<float> tempCalibUnc;
405  for (int crysID = 0; crysID < 8736; crysID++) {
406  tempCalib.push_back(CalibVsCrysID->GetBinContent(crysID + 1));
407  tempCalibUnc.push_back(CalibVsCrysID->GetBinError(crysID + 1));
408  }
409  ECLCrystalCalib* MuMuECalib = new ECLCrystalCalib();
410  MuMuECalib->setCalibVector(tempCalib, tempCalibUnc);
411  saveCalibration(MuMuECalib, "ECLCrystalEnergyMuMu");
412  B2RESULT("eclMuMuEAlgorithm: successfully stored ECLCrystalEnergyMuMu calibration constants");
413  }
414  }
415 
419  PeakVsCrysID->Write();
420  effSigVsCrysID->Write();
421  etaVsCrysID->Write();
422  constVsCrysID->Write();
423  normVsCrysID->Write();
424  fitLimitVsCrysID->Write();
425  StatusVsCrysID->Write();
426  hPeak->Write();
427  fracPeakUnc->Write();
428  nIterations->Write();
429  hStatus->Write();
430 
432  if (findExpValues) {
433  ExpEnergyperCrys->Write();
434  } else {
435  CalibVsCrysID->Write();
436  }
437  histfile->Close();
438 
441  dummy = (TH1F*)gROOT->FindObject("PeakVsCrysID"); delete dummy;
442  dummy = (TH1F*)gROOT->FindObject("effSigVsCrysID"); delete dummy;
443  dummy = (TH1F*)gROOT->FindObject("etaVsCrysID"); delete dummy;
444  dummy = (TH1F*)gROOT->FindObject("constVsCrysID"); delete dummy;
445  dummy = (TH1F*)gROOT->FindObject("normVsCrysID"); delete dummy;
446  dummy = (TH1F*)gROOT->FindObject("fitLimitVsCrysID"); delete dummy;
447  dummy = (TH1F*)gROOT->FindObject("StatusVsCrysID"); delete dummy;
448  dummy = (TH1F*)gROOT->FindObject("fitProbSame"); delete dummy;
449  dummy = (TH1F*)gROOT->FindObject("fracPeakUnc"); delete dummy;
450  dummy = (TH1F*)gROOT->FindObject("nIterations"); delete dummy;
451  dummy = (TH1F*)gROOT->FindObject("hStatus"); delete dummy;
452  dummy = (TH1F*)gROOT->FindObject("ExpEnergyperCrys"); delete dummy;
453  dummy = (TH1F*)gROOT->FindObject("CalibVsCrysID"); delete dummy;
454 
455 
458  if (storeConst == -1) {
459  B2RESULT("eclMuMuEAlgorithm performed fits but was not asked to store contants");
460  return c_Failure;
461  } else if (!DBsuccess) {
462  if (findExpValues) { B2RESULT("eclMuMuEAlgorithm: failed to store expected values"); }
463  else { B2RESULT("eclMuMuEAlgorithm: failed to store calibration constants"); }
464  return c_Failure;
465  }
466  return c_OK;
467 }
Belle2::ECL::eclMuMuEAlgorithm::notFit
int notFit
no fit performed
Definition: eclMuMuEAlgorithm.h:73
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::eclMuMuEAlgorithm::storeConst
int storeConst
controls which values are written to the database.
Definition: eclMuMuEAlgorithm.h:57
Belle2::ECL::eclMuMuEAlgorithm::tRatioMin
double tRatioMin
Fit range is adjusted so that fit at upper endpoint is between tRatioMin and tRatioMax of peak.
Definition: eclMuMuEAlgorithm.h:53
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::eclMuMuEAlgorithm::atLimit
int atLimit
a parameter is at the limit; fit not useable
Definition: eclMuMuEAlgorithm.h:70
Belle2::ECL::eclMuMuEAlgorithm::poorFit
int poorFit
low chi square; fit not useable
Definition: eclMuMuEAlgorithm.h:71
Belle2::ECL::eclMuMuEAlgorithm::findExpValues
bool findExpValues
if true, fits are used to find expected energy deposit for each crystal instead of the calibration co...
Definition: eclMuMuEAlgorithm.h:56
Belle2::ECL::eclMuMuEAlgorithm::cellIDHi
int cellIDHi
Last cellID to be fit.
Definition: eclMuMuEAlgorithm.h:50
Belle2::ECL::eclMuMuEAlgorithm::performFits
bool performFits
if false, input histograms are copied to output, but no fits are done.
Definition: eclMuMuEAlgorithm.h:55
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::ECL::eclMuMuEAlgorithm::calibrate
virtual EResult calibrate() override
..Run algorithm on events
Definition: eclMuMuEAlgorithm.cc:51
Belle2::ECL::eclMuMuEAlgorithm::cellIDLo
int cellIDLo
..Parameters to control Novosibirsk fit to energy deposited in each crystal by mu+mu- events
Definition: eclMuMuEAlgorithm.h:49
Belle2::CalibrationAlgorithm::c_Failure
@ c_Failure
Failed =3 in Python.
Definition: CalibrationAlgorithm.h:54
Belle2::ECL::eclMuMuEAlgorithm::noPeak
int noPeak
Novosibirsk component of fit is negligible; fit not useable.
Definition: eclMuMuEAlgorithm.h:72
Belle2::ECL::eclMuMuEAlgorithm::minEntries
int minEntries
All crystals to be fit must have at least minEntries events in the fit range.
Definition: eclMuMuEAlgorithm.h:51
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
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::ECL::eclMuMuEAlgorithm::fitOK
int fitOK
fit is OK
Definition: eclMuMuEAlgorithm.h:68
Belle2::ECL::eclMuMuEAlgorithm::tRatioMax
double tRatioMax
Fit range is adjusted so that fit at upper endpoint is between tRatioMin and tRatioMax of peak.
Definition: eclMuMuEAlgorithm.h:54
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47
Belle2::ECL::eclMuMuEAlgorithm::iterations
int iterations
fit reached max number of iterations, but is useable
Definition: eclMuMuEAlgorithm.h:69
Belle2::ECL::eclMuMuEAlgorithm::maxIterations
int maxIterations
no more than maxIteration iterations
Definition: eclMuMuEAlgorithm.h:52