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>
19 #include "TDirectory.h"
26 using namespace Calibration;
61 std::vector<double> eclLeakageFitParameters(TH1F* h,
const double& target)
65 double maxIntegral = h->Integral();
66 const int nBins = h->GetNbinsX();
72 std::vector<double> intVector;
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);
80 for (
int iLo = 2; iLo <= nBins; iLo++) {
81 for (
int iHi = iLo; iHi <= nBins; iHi++) {
84 double integral = intVector[iHi - 1] - intVector[iLo - 2];
87 if ((integral > target and (iHi - iLo) < (minHi - minLo)) or
88 (integral > target and (iHi - iLo) == (minHi - minLo) and integral > maxIntegral)
92 maxIntegral = integral;
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;
103 const int nfitBins = (1 + minHi - minLo);
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);
118 double eclLeakageNovo(Double_t* x, Double_t* par)
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;
128 if (TMath::Abs(tail) < 1.e-7)
129 qc = 0.5 * TMath::Power(((x[0] - peak) / width), 2);
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;
137 qc = 0.5 * (TMath::Power((log(qy) / tail), 2) + tail * tail);
142 return par[0] * exp(-qc);
149 int eclLeakageFitQuality(
const double& fitLo,
const double& fitHi,
const double& fitPeak,
const double& fitSigma,
150 const double& fitdEta,
const double& fitProb)
152 const double tolerance = 0.02;
153 const double redoFitProb = 1e-5;
154 const double minFitProb = 1e-9;
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;}
171 "Generate payload ECLLeakageCorrection from single photon MC recorded by eclLeakageCollectorModule"
184 B2INFO(
"starting eclLeakageAlgorithm");
188 auto inputParameters = getObjectPtr<TH1F>(
"inputParameters");
189 int lastBin = inputParameters->GetNbinsX();
190 double scale = inputParameters->GetBinContent(lastBin);
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.);
198 TFile* histfile =
new TFile(
"eclLeakageAlgorithm.root",
"recreate");
199 inputParameters->Write();
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.;
218 float generatedE[nLeakReg][nEnergies];
220 for (
int ireg = 0; ireg < nLeakReg; ireg++) {
221 B2INFO(
"Generated energies for ireg = " << ireg);
222 for (
int ie = 0; ie < nEnergies; ie++) {
224 generatedE[ireg][ie] = inputParameters->GetBinContent(bin);
225 B2INFO(
" " << ie <<
" " << generatedE[ireg][ie] <<
" GeV");
232 int iEnergiesMeV[nEnergies][nThetaID];
233 for (
int thID = 0; thID < nThetaID; thID++) {
235 if (thID >= firstBarrelID and thID <= lastBarrelID) {
237 }
else if (thID > lastBarrelID) {
240 for (
int ie = 0; ie < nEnergies; ie++) {
241 iEnergiesMeV[ie][thID] = (int)(1000.*generatedE[ireg][ie] + 0.5);
247 const double eFracLo = 0.4;
248 const double eFracHi = 1.5;
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++) {
255 double res_squared = 0.0001 + 0.064 / iEnergiesMeV[ie][thID];
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]);
267 TString statusString[7] = {
"good",
"refit",
"lowStat",
"lowProb",
"peakAtLimit",
"sigmaAtLimit",
"etaAtLimit"};
268 const double fracEnt[2] = {0.683, 0.5};
269 const double minEntries = 100.;
270 const double minMaxBin = 50.;
271 const double highMaxBin = 300.;
273 TF1* func =
new TF1(
"eclLeakageNovo", eclLeakageNovo, eFracLo, eFracHi, 4);
274 func->SetParNames(
"normalization",
"peak",
"effSigma",
"eta");
275 func->SetLineColor(kRed);
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);
287 tree->SetBranchAddress(
"nCrys", &
t_nCrys);
292 const int treeEntries = tree->GetEntries();
293 B2INFO(
"eclLeakageAlgorithm entries = " << treeEntries);
301 const int iRun = exprun[0].second;
302 const int iExp = exprun[0].first;
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");
322 existingThetaCorrection.Write();
323 existingPhiCorrection.Write();
324 existingCrysCorrection.Write();
332 TH1F* locError[nEnergies][nThetaID];
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.);
345 for (
int i = 0; i < treeEntries; i++) {
347 if (
t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0) {
354 const double minLocEntries = 49.5;
355 const double peakFrac = 0.025;
356 const double startOfPed = 10.0001;
358 float maxLocCut[nEnergies][nThetaID];
361 TH1F* locCutSummary =
new TH1F(
"locCutSummary",
"location cut for each thetaID/energy; xBin = thetaID + ie*nTheta", nbinX, 0,
363 TH1F* locCutEfficiency =
new TH1F(
"locCutEfficiency",
364 "Fraction of events passing location cut for each thetaID/energy; xBin = thetaID + ie*nTheta", nbinX, 0, nbinX);
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;}
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();
379 const int ix = 1 + thID + nThetaID * ie;
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);
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);
391 const double entries = locError[ie][thID]->GetEntries();
392 const double pass = locError[ie][thID]->Integral(1, ib - 1);
393 locCutEfficiency->SetBinContent(ix, pass / entries);
394 locCutEfficiency->SetBinError(ix, 0);
399 B2DEBUG(25,
" thID " << thID <<
" E " << iEnergiesMeV[ie][thID] <<
" ped " << pedestal <<
" threshold " << threshold <<
400 " location cut " << maxLocCut[ie][thID]);
403 title = locError[ie][thID]->GetTitle();
405 const double locEff = locCutEfficiency->GetBinContent(ix);
406 sCut.Form(
" cut %0.1f cm, eff %0.2f", maxLocCut[ie][thID], locEff);
408 locError[ie][thID]->SetTitle(title);
410 locError[ie][thID]->Write();
416 locCutSummary->Write();
417 locCutEfficiency->Write();
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";
434 hELabUncorr[ie][thID] =
new TH1F(name, title, 2 * nEfracBins[ie][thID], eFracLo, eFracHi);
440 for (
int i = 0; i < treeEntries; i++) {
444 bool goodReco =
t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0
446 if (not goodReco) {
continue;}
455 float peakUncorr[nEnergies][nThetaID];
456 float etaUncorr[nEnergies][nThetaID];
457 std::vector<TString> failedELabUncorr;
458 std::vector<int> statusELabUncorr;
459 int payloadStatus = 0;
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",
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];
476 double entries = hEnergy->Integral();
478 double genE = iEnergiesMeV[ie][thID] / 1000.;
479 bool fitHist = entries > minEntries;
485 double norm = hEnergy->GetMaximum();
486 double target = fracEnt[nIter] * entries;
487 std::vector<double> startingParameters;
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);
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();
505 double dEta = min((etaMax - eta), (eta - etaMin));
506 fitStatus = eclLeakageFitQuality(startingParameters[3], startingParameters[4], peak, effSigma, dEta, prob);
509 if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
519 statusOfhELabUncorr->Fill(fitStatus + 0.000001);
520 if (fitStatus == 2 or fitStatus >= 4) {
521 statusELabUncorr.push_back(fitStatus);
522 failedELabUncorr.push_back(hEnergy->GetName());
529 peakUncorr[ie][thID] = peak;
530 etaUncorr[ie][thID] = eta;
533 if (entries > minEntries) {
535 hELabUncorr[ie][thID]->Write();
542 statusOfhELabUncorr->Write();
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());
552 if (payloadStatus != 0) {
553 B2ERROR(
"ecLeakageAlgorithm: fit to hELabUncorr failed. ");
560 for (
int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
561 for (
int ie = 0; ie < nEnergies; ie++) {
562 if (peakUncorr[ie][thID] < 0.) {
564 for (
int nextID = thID - 1; thID >= firstUsefulThID; thID--) {
565 if (peakUncorr[ie][nextID] > 0.) {
566 peakUncorr[ie][thID] = peakUncorr[ie][nextID];
571 for (
int nextID = thID + 1; thID <= lastUsefulThID; thID++) {
572 if (peakUncorr[ie][nextID] > 0.) {
573 peakUncorr[ie][thID] = peakUncorr[ie][nextID];
581 if (peakUncorr[ie][thID] < 0.) {
582 B2ERROR(
"ecLeakageAlgorithm: unable to correct hELabUncorr for thetaID " << thID <<
" ie " << ie);
599 TString dirName[nDir] = {
"theta",
"phiMech",
"phiNoMech"};
601 TH1F* eFracPosition[nEnergies][nThetaID][nDir][nPositions];
602 std::vector<TString> failedeFracPosition;
603 std::vector<int> statuseFracPosition;
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 +
616 eFracPosition[ie][thID][idir][ipos] =
new TH1F(name, title, nEfracBins[ie][thID], eFracLo, eFracHi);
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,
626 TH1F* probOfeFracPosition =
new TH1F(
"probOfeFracPosition",
"fit probability of eFrac fits for each position;probability", 100, 0,
628 TH1F* maxOfeFracPosition =
new TH1F(
"maxOfeFracPosition",
"max entries of eFrac histograms;maximum bin content", 100, 0, 1000);
632 for (
int i = 0; i < treeEntries; i++) {
636 bool goodReco =
t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0
638 if (not goodReco) {
continue;}
651 float positionCorrection[nEnergies][nThetaID][nDir][nPositions];
652 float positionCorrectionUnc[nEnergies][nThetaID][nDir][nPositions];
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);
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++;}
670 double correction = sqrt(peakUncorr[ie][thID]);
671 double corrUnc = 0.05;
675 if (hEnergy->GetEntries() < 0.5) {fitStatus = -2;}
677 double entries = hEnergy->Integral();
685 double target = fracEnt[nIter] * entries;
686 std::vector<double> startingParameters;
687 bool findParm =
true;
689 startingParameters = eclLeakageFitParameters(hEnergy, target);
692 if (hEnergy->GetMaximum()<minMaxBin and startingParameters[5]>11.9) {
698 double norm = hEnergy->GetMaximum();
699 double fitLow = startingParameters[3];
700 double fitHigh = startingParameters[4];
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);
709 if (hEnergy->GetMaximum() < highMaxBin) {
710 func->SetParLimits(3, etaFix, etaFix);
712 func->SetParLimits(3, etaMin, etaMax);
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();
726 double dEta = min((etaMax - eta), (eta - etaMin));
727 fitStatus = eclLeakageFitQuality(fitLow, fitHigh, peak, effSigma, dEta, prob);
730 if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
735 }
else if (fitStatus <= 3) {
739 correction = peak / sqrt(peakUncorr[ie][thID]);
740 corrUnc = func->GetParError(1) / sqrt(peakUncorr[ie][thID]);
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());
752 if (prob > 0 and fitStatus <= 3) {
754 thetaCorSummary->Fill(correction);
756 phiCorSummary->Fill(correction);
761 if (fitStatus >= 2) {
762 statuseFracPosition.push_back(fitStatus);
763 failedeFracPosition.push_back(hEnergy->GetName());
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());
784 statusOfeFracPosition->Write();
785 probOfeFracPosition->Write();
786 maxOfeFracPosition->Write();
787 thetaCorSummary->Write();
788 phiCorSummary->Write();
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";
805 hEPos[ie][thID] =
new TH1F(name, title, 2 * nEfracBins[ie][thID], eFracLo, eFracHi);
811 for (
int i = 0; i < treeEntries; i++) {
815 bool goodReco =
t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0
817 if (not goodReco) {
continue;}
834 float etaEpos[nEnergies][nThetaID];
835 std::vector<TString> failedEPos;
836 std::vector<int> statusEPos;
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",
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];
852 double entries = hEnergy->Integral();
860 double norm = hEnergy->GetMaximum();
861 double target = fracEnt[nIter] * entries;
862 std::vector<double> startingParameters;
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);
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();
880 double dEta = min((etaMax - eta), (eta - etaMin));
881 fitStatus = eclLeakageFitQuality(startingParameters[3], startingParameters[4], peak, effSigma, dEta, prob);
884 if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
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; }
902 etaEpos[ie][thID] = eta;
905 if (entries > minEntries) {
907 hEPos[ie][thID]->Write();
914 statusOfhEPos->Write();
916 if (payloadStatus != 0) {
917 B2ERROR(
"ecLeakageAlgorithm: fit to hEPos failed. ");
931 TH1F* ePosnCry[nEnergies][nThetaID][maxN + 1];
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);
947 for (
int i = 0; i < treeEntries; i++) {
951 bool goodReco =
t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0
953 if (not goodReco) {
continue;}
968 float nCrysCorrection[nEnergies][nThetaID][maxN + 1];
969 float nCrysCorrectionUnc[nEnergies][nThetaID][maxN + 1];
973 int maxNCry[nEnergies][nThetaID];
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,
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++) {
987 TH1F* hEnergy = (TH1F*)ePosnCry[ie][thID][in];
988 if (hEnergy->Integral() > minEntries) {nHistToFit++;}
991 double correction = -1.;
992 double corrUnc = 0.05;
994 double entries = hEnergy->Integral();
1004 if (entries > maxIntegral) {
1005 maxIntegral = entries;
1006 maxNCry[ie][thID] = in;
1014 double target = fracEnt[nIter] * entries;
1015 std::vector<double> startingParameters;
1016 bool findParm =
true;
1018 startingParameters = eclLeakageFitParameters(hEnergy, target);
1021 if (hEnergy->GetMaximum()<minMaxBin and startingParameters[5]>11.9) {
1027 double norm = hEnergy->GetMaximum();
1028 double fitLow = startingParameters[3];
1029 double fitHigh = startingParameters[4];
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);
1038 if (hEnergy->GetMaximum() < highMaxBin) {
1039 func->SetParLimits(3, etaFix, etaFix);
1041 func->SetParLimits(3, etaMin, etaMax);
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();
1055 double dEta = min((etaMax - eta), (eta - etaMin));
1056 fitStatus = eclLeakageFitQuality(fitLow, fitHigh, peak, effSigma, dEta, prob);
1059 if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
1064 }
else if (fitStatus <= 3) {
1066 corrUnc = func->GetParError(1);
1071 statusOfePosnCry->Fill(fitStatus + 0.00001);
1072 nCrysCorrection[ie][thID][in] = correction;
1073 nCrysCorrectionUnc[ie][thID][in] = corrUnc;
1078 int highestN = maxNCry[ie][thID];
1079 ePosnCry[ie][thID][highestN]->Write();
1085 statusOfePosnCry->Write();
1091 for (
int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
1092 for (
int ie = 0; ie < nEnergies; ie++) {
1096 int highestN = maxNCry[ie][thID];
1097 if (nCrysCorrection[ie][thID][highestN] <= 0) {
1098 B2ERROR(
"eclLeakageAlgorithm: no nCrys correction for ie " << ie <<
" thetaID " << thID);
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];
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];
1130 for (
int thID = firstUsefulThID - 1; thID >= 0; thID--) {
1131 for (
int ie = 0; ie < nEnergies; ie++) {
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];
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];
1150 for (
int thID = lastUsefulThID + 1; thID < nThetaID; thID++) {
1151 for (
int ie = 0; ie < nEnergies; ie++) {
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];
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];
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]));
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];
1192 TH2F thetaCorrection(
"thetaCorrection",
"Theta correction vs bin;bin = thetaID + 69*energyBin;local theta bin", nbinX, 0, nbinX,
1193 nPositions, 0, nPositions);
1195 for (
int ie = 0; ie < nEnergies; ie++) {
1196 for (
int thID = 0; thID < nThetaID; thID++) {
1198 for (
int ipos = 0; ipos < nPositions; ipos++) {
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);
1210 TH2F phiCorrection(
"phiCorrection",
"Phi correction vs bin;bin = thetaID + 69*energyBin;local phi bin", 2 * nbinX, 0, 2 * nbinX,
1211 nPositions, 0, nPositions);
1213 for (
int ie = 0; ie < nEnergies; ie++) {
1214 for (
int thID = 0; thID < nThetaID; thID++) {
1216 for (
int ipos = 0; ipos < nPositions; ipos++) {
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);
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);
1236 TH2F nCrystalCorrection(
"nCrystalCorrection",
"nCrys correction vs bin;bin = thetaID + 69*energyBin;nCrys", nbinX, 0, nbinX,
1237 maxN + 1, 0, maxN + 1);
1240 for (
int ie = 0; ie < nEnergies; ie++) {
1241 for (
int thID = 0; thID < nThetaID; thID++) {
1243 for (
int in = 0; in <= maxN; in++) {
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);
1260 thetaCorrection.Write();
1261 phiCorrection.Write();
1262 nCrystalCorrection.Write();
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];
1273 int thIDReg[nLeakReg];
1278 for (
int ireg = 0; ireg < nLeakReg; ireg++) {
1279 for (
int ie = 0; ie < nEnergies; ie++) {
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]];
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);
1297 for (
int i = 0; i < treeEntries; i++) {
1301 bool goodReco =
t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0
1303 if (not goodReco) {
continue;}
1314 float posCor = thetaPosCor * phiPosCor;
1318 float corrTrue = posCor * nCrysCor;
1323 float corrNonCrys = 0.96;
1324 for (
int iter = 0; iter < 2; iter++) {
1329 float logEnergy = log(energyRaw / corrNonCrys);
1331 if (logEnergy < leakLogE[
t_region][0]) {
1333 }
else if (logEnergy > leakLogE[
t_region][nEnergies - 1]) {
1334 ie0 = nEnergies - 2;
1336 while (logEnergy > leakLogE[
t_region][ie0 + 1]) {ie0++;}
1345 corrNonCrys = cor0 + (cor1 - cor0) * (logEnergy - leakLogE[
t_region][ie0]) / (leakLogE[
t_region][ie0 + 1] -
1352 float corrMeasured = 0.96;
1353 for (
int iter = 0; iter < 2; iter++) {
1358 float logEnergy = log(energyRaw / corrMeasured);
1360 if (logEnergy < leakLogE[
t_region][0]) {
1362 }
else if (logEnergy > leakLogE[
t_region][nEnergies - 1]) {
1363 ie0 = nEnergies - 2;
1365 while (logEnergy > leakLogE[
t_region][ie0 + 1]) {ie0++;}
1375 corrMeasured = cor0 + (cor1 - cor0) * (logEnergy - leakLogE[
t_region][ie0]) / (leakLogE[
t_region][ie0 + 1] -
1392 float peakEnergy[nLeakReg][nEnergies][nResType];
1393 float energyRes[nLeakReg][nEnergies][nResType];
1396 for (
int ireg = 0; ireg < nLeakReg; ireg++) {
1397 for (
int ie = 0; ie < nEnergies; ie++) {
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];
1406 bool fitHist = entries > minEntries;
1412 double norm = hEnergy->GetMaximum();
1413 double target = fracEnt[nIter] * entries;
1414 std::vector<double> startingParameters;
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);
1425 int minLo = hEnergy->GetXaxis()->FindBin(startingParameters[3] + 0.0001);
1426 int minHi = hEnergy->GetXaxis()->FindBin(startingParameters[4] - 0.0001);
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);
1435 res68 = 0.5 * dx * (1 + minHi - minLo - max(subLo, subHi));
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();
1449 double dEta = min((etaMax - eta), (eta - etaMin));
1450 int fitStatus = eclLeakageFitQuality(startingParameters[3], startingParameters[4], peak, effSigma, dEta, prob);
1453 if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
1462 peakEnergy[ireg][ie][ires] = peak;
1463 energyRes[ireg][ie][ires] = res68;
1474 int nresBins = nEnergies * nLeakReg * (nResType + 1);
1475 TH1F* peakSummary =
new TH1F(
"peakSummary",
"Peak E/Etrue for each method, region, energy;Energy energy point", nresBins, 0,
1477 TH1F* resolutionSummary =
new TH1F(
"resolutionSummary",
"Resolution/peak for each method, region, energy;Energy energy point",
1478 nresBins, 0, nEnergies);
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]);}
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++) {
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.);
1497 B2INFO(
" " << ires <<
" " << energyRes[ireg][ie][ires] / peakEnergy[ireg][ie][ires]);
1505 peakSummary->Write();
1506 resolutionSummary->Write();
1528 B2INFO(
"eclLeakageAlgorithm: successfully stored payload ECLLeakageCorrections");
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.
Singleton class to cache database objects.
static DataStore & Instance()
Instance of singleton Store.
void setInitializeActive(bool active)
Setter for m_initializeActive.
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_thetaID
thetaID of photon
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.
bool construct(Args &&... params)
Construct an object of type T in this StoreObjPtr, using the provided constructor arguments.
static DBStore & Instance()
Instance of a singleton DBStore.
void update()
Updates all objects that are outside their interval of validity.
Abstract base class for different kinds of events.