Belle II Software development
eclNOptimalAlgorithm.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
9/* Own header. */
10#include <ecl/calibration/eclNOptimalAlgorithm.h>
11
12/* ECL headers. */
13#include <ecl/calibration/tools.h>
14#include <ecl/dataobjects/ECLElementNumbers.h>
15#include <ecl/dbobjects/ECLnOptimal.h>
16
17/* Basf2 headers. */
18#include <framework/database/DBObjPtr.h>
19#include <framework/database/DBStore.h>
20#include <framework/dataobjects/EventMetaData.h>
21#include <framework/datastore/DataStore.h>
22#include <framework/datastore/StoreObjPtr.h>
23
24/* ROOT headers. */
25#include <TF1.h>
26#include <TH2F.h>
27#include <TMath.h>
28#include <TROOT.h>
29
30using namespace std;
31using namespace Belle2;
32using namespace ECL;
33using namespace Calibration;
34
35//-----------------------------------------------------------------------------------
36//-----------------------------------------------------------------------------------
37//..Structure
38// 1. Setup
39// - from collector histograms, read in job parameters and group number for each crystal
40// - get the existing payload, which we will use as the starting point for nOptimal
41// - write out histograms from collector and from the existing payload
42//
43// 2. Create histograms to store many items from many fits
44//
45// 3. Loop over each energy and group to find the optimal number of crystals
46// - start by checking the current value of nOptimal, then check smaller and larger
47// values of number of summed crystals. For value, find:
48// - peak contained energy divided by generated energy by fitting the histogram
49// range that includes 50% of events;
50// - corresponding bias = sum of ECLCalDigits minus mc true energy;
51// - resolution = one-half of minimum range that includes 68.3% of events
52// - after finding the value of nOptimal with best resolution, also find the
53// resolution for the current reconstruction.
54// - store everything in diagnostic histograms.
55//
56// 4. Having found nOptimal for that group/energy, find the bias and peak fractional
57// contained energy for adjacent energies. This enables energy interpolation
58// when the payloads are used.
59//
60// 5. Generate and store the actual payload
61
62
63//-----------------------------------------------------------------------------------
64//-----------------------------------------------------------------------------------
65//..Some helper functions
66
67//-----------------------------------------------------------------------------------
68//..Fit the minimum range of bins that contain the specified number of events.
69std::vector<int> eclNOptimalFitRange(TH1D* h, const double& fraction)
70{
71 const double target = fraction * h->GetEntries();
72 const int nBins = h->GetNbinsX();
73
74 //..Start at the histogram maximum
75 int iLo = h->GetMaximumBin();
76 int iHi = iLo;
77 double sum = h->Integral(iLo, iHi);
78
79 //..Add one bin at a time
80 while (sum < target and (iLo > 1 or iHi < nBins)) {
81 double nextLo = 0.;
82 if (iLo > 1) {nextLo = h->GetBinContent(iLo - 1);}
83 double nextHi = 0.;
84 if (iHi < nBins) {nextHi = h->GetBinContent(iHi + 1);}
85 if (nextLo > nextHi) {
86 sum += nextLo;
87 iLo--;
88 } else {
89 sum += nextHi;
90 iHi++;
91 }
92 }
93
94 std::vector<int> iBins;
95 iBins.push_back(iLo);
96 iBins.push_back(iHi);
97 return iBins;
98
99}
100
101//-----------------------------------------------------------------------------------
102//..Resolution is the minimum range that contains 68.3% of entries
103double eclNOptimalResolution(TH1D* h, int& iLo75, int& iHi75)
104{
105
106 //..Search among the bin range that contains 75% of events for the smallest
107 // range that contains at least 68.3%. If more than one, pick the one with
108 // the most events.
109 const int nBins = h->GetNbinsX();
110 const double target = 0.683 * h->GetEntries();
111
112 //..Copy the histogram contents for speed of access
113 // Recall that the first histogram bin is 1, not 0.
114 std::vector<double> intVector;
115 intVector.push_back(h->GetBinContent(1));
116 for (int iL = 2; iL <= nBins; iL++) {
117 double nextIntegral = intVector[iL - 2] + h->GetBinContent(iL);
118 intVector.push_back(nextIntegral);
119 }
120
121 //..Now search all possible ranges
122 int iLo = iLo75;
123 int iHi = iHi75;
124 double maxIntegral = intVector[iHi - 1] - intVector[iLo - 2];
125 for (int iL = iLo75; iL <= iHi75; iL++) {
126 for (int iH = iL; iH <= iHi75; iH++) {
127
128 // sum[iL, iH] = sum[1, iH] - sum[1,iL-1] = intVector[iH-1] - intVector[iL-2]
129 double integral = intVector[iH - 1] - intVector[iL - 2];
130 if ((integral > target and (iH - iL) < (iHi - iLo)) or
131 (integral > target and (iH - iL) == (iHi - iLo) and integral > maxIntegral)
132 ) {
133 iLo = iL;
134 iHi = iH;
135 maxIntegral = integral;
136 }
137 }
138 }
139
140 //..Refine the range to exactly 68.3% by subtracting a fraction of the first or last bin.
141 // Use fit to estimate distribution of events across the bin.
142 const double overage = h->Integral(iLo, iHi) - target;
143 const double dx = (h->GetXaxis()->GetXmax() - h->GetXaxis()->GetXmin()) / nBins;
144
145 //..Fraction of the first bin we could exclude
146 double xLow = h->GetBinLowEdge(iLo);
147 double fracEntriesToExcludeLo = overage / h->GetBinContent(iLo);
148 double fracBinToExcludeLo = fracEntriesToExcludeLo;
149
150 //..Use the slope from the fit extrapolated to this point unless it is 0
151 TF1* func = (TF1*)h->GetFunction("eclNOptimalNovo");
152 double f0 = func->Eval(xLow);
153 double f1 = func->Eval(xLow + dx);
154 if (abs(f1 - f0) > 1.) {
155 fracBinToExcludeLo = (sqrt(f0 * f0 + fracEntriesToExcludeLo * (f1 * f1 - f0 * f0)) - f0) / (f1 - f0);
156 }
157
158 //..Last bin. In this case, we want the integral from the low edge of the bin to
159 // be the fraction of events to keep in the range, not the fraction to exclude.
160 double xHigh = h->GetBinLowEdge(iHi + 1);
161 f0 = func->Eval(xHigh - dx);
162 f1 = func->Eval(xHigh);
163 double fracEntriesToExcludeHi = overage / h->GetBinContent(iHi);
164 double fracEntriesToInclude = 1. - fracEntriesToExcludeHi;
165 double fracBinToInclude = fracEntriesToInclude;
166 if (abs(f1 - f0) > 1.) {
167 fracBinToInclude = (sqrt(f0 * f0 + fracEntriesToInclude * (f1 * f1 - f0 * f0)) - f0) / (f1 - f0);
168 }
169 double fracBinToExcludeHi = 1. - fracBinToInclude;
170
171 //..Exclude the first bin or last bin fraction, whichever gives better resolution
172 double rangeBins = 1 + iHi - iLo - max(fracBinToExcludeLo, fracBinToExcludeHi);
173 return 0.5 * dx * rangeBins;
174
175}
176
177//--------------------------------------------------------------------
178//..Novosibirsk; H. Ikeda et al. / NIM A 441 (2000) 401-426
179double eclNOptimalNovo(const double* x, const double* par)
180{
181
182 double peak = par[1];
183 double width = par[2];
184 double sln4 = sqrt(log(4));
185 double y = par[3] * sln4;
186 double tail = -log(y + sqrt(1 + y * y)) / sln4;
187 double qc = 0.;
188
189 if (TMath::Abs(tail) < 1.e-7)
190 qc = 0.5 * TMath::Power(((x[0] - peak) / width), 2);
191 else {
192 double qa = tail * sqrt(log(4.));
193 double qb = sinh(qa) / qa;
194 double qx = (x[0] - peak) / width * qb;
195 double qy = 1. + tail * qx;
196
197 if (qy > 1.E-7)
198 qc = 0.5 * (TMath::Power((log(qy) / tail), 2) + tail * tail);
199 else
200 qc = 15.0;
201 }
202
203 return par[0] * exp(-qc);
204}
205
206//-----------------------------------------------------------------------------------
207//-------------------------------------------------------------------------------
208//..Constructor
209
211{
213 "Use single gamma MC to find the optimal number of crystals to sum for cluster energy"
214 );
215}
216
217//-------------------------------------------------------------------------------
218//..Calibrate
220{
221
222 //-----------------------------------------------------------------------------------
223 //-----------------------------------------------------------------------------------
224 //..Algorithm set up
225
226 //..Get started
227 B2INFO("Starting eclNOptimalAlgorithm");
228 gROOT->SetBatch();
229
230 //-----------------------------------------------------------------------------------
231 //..Read in parameter histograms created by the collector and fix normalization
232 auto inputParameters = getObjectPtr<TH1F>("inputParameters");
233 const int lastBin = inputParameters->GetNbinsX();
234 const double scale = inputParameters->GetBinContent(lastBin); // number of times inputParameters was filled
235 for (int ib = 1; ib < lastBin; ib++) {
236 double param = inputParameters->GetBinContent(ib);
237 inputParameters->SetBinContent(ib, param / scale);
238 inputParameters->SetBinError(ib, 0.);
239 }
240
241 //..Also read in and normalize group number for each cellID
242 auto groupNumberOfEachCellID = getObjectPtr<TH1F>("groupNumberOfEachCellID");
243 for (int ic = 1; ic <= ECLElementNumbers::c_NCrystals; ic++) {
244 const double groupNum = groupNumberOfEachCellID->GetBinContent(ic);
245 groupNumberOfEachCellID->SetBinContent(ic, groupNum / scale);
246 groupNumberOfEachCellID->SetBinError(ic, 0.);
247 }
248
249 //..Couple of diagnostic histograms
250 auto entriesPerThetaIdEnergy = getObjectPtr<TH2F>("entriesPerThetaIdEnergy");
251 auto mcEnergyDiff = getObjectPtr<TH2F>("mcEnergyDiff");
252
253 //..Write these to disk.
254 TFile* histFile = new TFile("eclNOptimalAlgorithm.root", "recreate");
255 inputParameters->Write();
256 groupNumberOfEachCellID->Write();
257 entriesPerThetaIdEnergy->Write();
258 mcEnergyDiff->Write();
259
260 //-----------------------------------------------------------------------------------
261 //..Parameters from the inputParameters histogram
262 const int nCrystalGroups = (int)(inputParameters->GetBinContent(1) + 0.5);
263 const int nEnergies = (int)(inputParameters->GetBinContent(2) + 0.5);
264 const int nLeakReg = 3; // 3 regions forward barrel backward
265 auto generatedE = create2Dvector<float>(nLeakReg, nEnergies);
266 int bin = 2; // bin 1 = nCrystalGroups, bin 2 = nEnergies, bin 3 = first energy
267 for (int ireg = 0; ireg < nLeakReg; ireg++) {
268 for (int ie = 0; ie < nEnergies; ie++) {
269 bin++;
270 generatedE[ireg][ie] = inputParameters->GetBinContent(bin);
271 }
272 }
273 B2INFO("nCrystalGroups = " << nCrystalGroups);
274
275 //-----------------------------------------------------------------------------------
276 //..Experiment and run number, for reading existing payload from database
277 const auto exprun = getRunList();
278 const int iExp = exprun[0].first;
279 const int iRun = exprun[0].second;
280 B2INFO("Experiment, run = " << iExp << ", " << iRun);
281
283 // simulate the initialize() phase where we can register objects in the DataStore
285 evtPtr.registerInDataStore();
287 // now construct the event metadata
288 evtPtr.construct(1, iRun, iExp);
289 // and update the database contents
290 DBStore& dbstore = DBStore::Instance();
291 dbstore.update();
292 // this is only needed it the payload might be intra-run dependent,
293 // that is if it might change during one run as well
294 dbstore.updateEvent();
295
296 //-----------------------------------------------------------------------------------
297 //..Get existing payload
298 B2INFO("eclNOptimalAlgorithm: reading previous payload");
299 DBObjPtr<ECLnOptimal> existingECLnOptimal;
300 std::vector<float> eUpperBoundariesFwdPrev = existingECLnOptimal->getUpperBoundariesFwd();
301 std::vector<float> eUpperBoundariesBrlPrev = existingECLnOptimal->getUpperBoundariesBrl();
302 std::vector<float> eUpperBoundariesBwdPrev = existingECLnOptimal->getUpperBoundariesBwd();
303 TH2F nOptimal2DPrev = existingECLnOptimal->getNOptimal();
304 nOptimal2DPrev.SetName("nOptimal2DPrev");
305
306 const int nPrevE = eUpperBoundariesFwdPrev.size();
307 B2INFO("Energy upper boundaries from previous payload for each region");
308 for (int ie = 0; ie < nPrevE; ie++) {
309 B2INFO(" " << ie << " " << eUpperBoundariesFwdPrev[ie] << " " << eUpperBoundariesBrlPrev[ie] << " " << eUpperBoundariesBwdPrev[ie]);
310 }
311
312 //..Write out
313 histFile->cd();
314 nOptimal2DPrev.Write();
315
316 //-----------------------------------------------------------------------------------
317 //..Use existing payload to get the initial nOptimal value for each group
318 std::vector< std::vector<int> > initialNOptimal;
319 initialNOptimal.resize(nCrystalGroups, std::vector<int>(nEnergies, 0));
320
321 //..Couple of histograms of relevant quantities
322 TH2F* nInitialPerGroup = new TH2F("nInitialPerGroup", "initial nCrys, energy point vs group number;group;energy point",
323 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
324 TH2F* generatedEPerGroup = new TH2F("generatedEPerGroup", "Generated energy, energy point vs group number;group;energy point",
325 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
326 TH1F* firstCellIDPerGroup = new TH1F("firstCellIDPerGroup", "First cellID of group;group", nCrystalGroups, 0., nCrystalGroups);
327 TH1F* nCrystalsPerGroup = new TH1F("nCrystalsPerGroup", "Number of crystals per group;group", nCrystalGroups, 0., nCrystalGroups);
328 std::vector<int> crystalsPerGroup;
329 crystalsPerGroup.resize(nCrystalGroups, 0);
330
331 const int nCrystalsTotal = ECLElementNumbers::c_NCrystals;
332 for (int cellID = 1; cellID <= nCrystalsTotal; cellID++) {
333
334 //..Group number of this cellID
335 int iGroup = (int)(0.5 + groupNumberOfEachCellID->GetBinContent(cellID));
336
337 //..Energy boundaries of previous payload, which depend on the ECL region
338 std::vector<float> eUpperBoundariesPrev = eUpperBoundariesBrlPrev;
339 int iRegion = 1; // barrel
340 if (ECLElementNumbers::isForward(cellID)) {
341 eUpperBoundariesPrev = eUpperBoundariesFwdPrev;
342 iRegion = 0;
343 } else if (ECLElementNumbers::isBackward(cellID)) {
344 eUpperBoundariesPrev = eUpperBoundariesBwdPrev;
345 iRegion = 2;
346 }
347
348 //..For each test energy, get the nOptimal for this cellID from previous payload.
349 // iEnergy is the energy bin of the previous payload, which is equal to ie if
350 // the test energies have not changed.
351 for (int ie = 0; ie < nEnergies; ie++) {
352 float energy = generatedE[iRegion][ie];
353 int iEnergy = 0;
354 while (energy > eUpperBoundariesPrev[iEnergy]) {iEnergy++;}
355 initialNOptimal[iGroup][ie] = (int)(0.5 + nOptimal2DPrev.GetBinContent(iGroup + 1, iEnergy + 1));
356 if (ie != iEnergy) {
357 B2INFO("ie iEnergy mismatch: cellID " << cellID << " ie " << ie << " iEnergy " << iEnergy);
358 }
359
360 //..Store in diagnostic histograms
361 generatedEPerGroup->SetBinContent(iGroup + 1, ie + 1, energy);
362 nInitialPerGroup->SetBinContent(iGroup + 1, ie + 1, initialNOptimal[iGroup][ie]);
363 }
364
365 //..Count crystals in this group
366 if (crystalsPerGroup[iGroup] == 0) {
367 firstCellIDPerGroup->SetBinContent(iGroup + 1, cellID);
368 firstCellIDPerGroup->SetBinError(iGroup + 1, 0.);
369 }
370 crystalsPerGroup[iGroup]++;
371 }
372
373 //..Write out these histograms
374 for (int ig = 0; ig < nCrystalGroups; ig++) {
375 nCrystalsPerGroup->SetBinContent(ig + 1, crystalsPerGroup[ig]);
376 nCrystalsPerGroup->SetBinError(ig + 1, 0.);
377 }
378 histFile->cd();
379 nInitialPerGroup->Write();
380 generatedEPerGroup->Write();
381 firstCellIDPerGroup->Write();
382 nCrystalsPerGroup->Write();
383 B2INFO("Wrote initial diagnostic histograms");
384
385 //-----------------------------------------------------------------------------------
386 //..Write out all the 2D histograms
387 for (int ig = 0; ig < nCrystalGroups; ig++) {
388 for (int ie = 0; ie < nEnergies; ie++) {
389 std::string name = "eSum_" + std::to_string(ig) + "_" + std::to_string(ie);
390 auto eSum = getObjectPtr<TH2F>(name);
391 name = "biasSum_" + std::to_string(ig) + "_" + std::to_string(ie);
392 auto biasSum = getObjectPtr<TH2F>(name);
393 histFile->cd();
394 eSum->Write();
395 biasSum->Write();
396 }
397 }
398 B2INFO("Wrote eSum and biasSum histograms");
399
400 //-----------------------------------------------------------------------------------
401 //-----------------------------------------------------------------------------------
402 //..Find nOptimal and corresponding bias for each group/energy point
403
404 //..Histograms to store fit results
405 TH2F* nOptimalPerGroup = new TH2F("nOptimalPerGroup", "nOptimal;group;energy point", nCrystalGroups, 0., nCrystalGroups, nEnergies,
406 0., nEnergies);
407 TH2F* changeInNOptimal = new TH2F("changeInNOptimal", "nOptimal minus nInitial;group;energy point", nCrystalGroups, 0.,
408 nCrystalGroups, nEnergies, 0., nEnergies);
409 TH2F* binsInFit = new TH2F("binsInFit", "Number of bins used in Novo fit, energy point vs group number;group;energy point",
410 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
411 TH2F* maxEntriesPerHist = new TH2F("maxEntriesPerHist",
412 "Max entries in histogram bin, energy point vs group number;group;energy point", nCrystalGroups, 0., nCrystalGroups, nEnergies, 0.,
413 nEnergies);
414 TH2F* peakPerGroup = new TH2F("peakPerGroup", "peak energy, energy point vs group number;group;energy point", nCrystalGroups, 0.,
415 nCrystalGroups, nEnergies, 0., nEnergies);
416 TH2F* effSigmaPerGroup = new TH2F("effSigmaPerGroup", "fit effective sigma, energy point vs group number;group;energy point",
417 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
418 TH2F* etaPerGroup = new TH2F("etaPerGroup", "fit eta, energy point vs group number;group;energy point", nCrystalGroups, 0.,
419 nCrystalGroups, nEnergies, 0., nEnergies);
420 TH2F* fitProbPerGroup = new TH2F("fitProbPerGroup", "fit probability, energy point vs group number;group;energy point",
421 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
422 TH2F* resolutionPerGroup = new TH2F("resolutionPerGroup", "resolution/(peak-bias), energy point vs group number;group;energy point",
423 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
424 TH2F* fracContainedEPerGroup = new TH2F("fracContainedEPerGroup",
425 "peak fraction of energy contained in nOpt crystals;group;energy point", nCrystalGroups, 0., nCrystalGroups, nEnergies, 0.,
426 nEnergies);
427
428 TH2F* biasPerGroup = new TH2F("biasPerGroup", "bias (GeV), energy point vs group number;group;energy point", nCrystalGroups, 0.,
429 nCrystalGroups, nEnergies, 0., nEnergies);
430 TH2F* sigmaBiasPerGroup = new TH2F("sigmaBiasPerGroup", "sigma bias (GeV), energy point vs group number;group;energy point",
431 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
432
433 TH2F* existingResolutionPerGroup = new TH2F("existingResolutionPerGroup",
434 "existing resolution/peak, energy point vs group number;group;energy point", nCrystalGroups, 0., nCrystalGroups, nEnergies, 0.,
435 nEnergies);
436
437 //-----------------------------------------------------------------------------------
438 //..Loop over all groups and energy points
439
440 //..Fit function Novosibirsk (xmin, xmax, nparameters)
441 TF1* func = new TF1("eclNOptimalNovo", eclNOptimalNovo, 0.4, 1.4, 4);
442 func->SetLineColor(kRed);
443
444 for (int ig = 0; ig < nCrystalGroups; ig++) {
445 for (int ie = 0; ie < nEnergies; ie++) {
446
447 //..Read in the energy sum vs nCrys for the group and energy point
448 std::string name = "eSum_" + std::to_string(ig) + "_" + std::to_string(ie);
449 auto eSum = getObjectPtr<TH2F>(name);
450
451 //..Also the corresponding bias histogram
452 std::string biasName = "biasSum_" + std::to_string(ig) + "_" + std::to_string(ie);
453 auto biasSum = getObjectPtr<TH2F>(biasName);
454
455 //..And the generated energy
456 const double genEnergy = generatedEPerGroup->GetBinContent(ig + 1, ie + 1);
457
458 //-----------------------------------------------------------------------------------
459 //..Extract the eSum distributions for various nCrys (=nCrysSumToFit) and fit them
460 // to find the best resolution.
461 const int nCrysBins = eSum->GetNbinsX();
462 const int nCrysMax = nCrysBins - 1; // the last bin is eSum for existing reco
463
464 //..Items to be found and stored in the following while loop
465 TH1D* eSumOpt{nullptr}; // 1D projection with the best resolution
466 int nFitBins = 0; // keep track of how many bins are included in the fit
467 double bestFractionalResolution = 999.;
468 double existingFractionalResolution = 999.;
469 double biasForNOpt = 0.; // bias from backgrounds for optimal nCrys
470 double biasSigmaForNOpt = 0.; // sigma bias
471 int nOpt = 0;
472
473 //..Start with the previous optimal value, then try fewer, if possible, then
474 // try more, if possible. End by finding the resolution for the
475 // reconstruction used to produce the single photon MC samples.
476 const int initialnCrysSumToFit = initialNOptimal[ig][ie];
477 int nCrysSumToFit = initialnCrysSumToFit;
478 while (nCrysSumToFit > 0) {
479
480 TH1D* hEnergy = (TH1D*)eSum->ProjectionY("hEnergy", nCrysSumToFit, nCrysSumToFit);
481 TString newName = name + "_" + std::to_string(nCrysSumToFit);
482 hEnergy->SetName(newName);
483
484 //------------------------------------------------------------------------------
485 //..Check stats, and rebin as required
486 const double minESumEntries = 800.;
487 if (hEnergy->GetEntries() < minESumEntries) {
488 B2INFO("Insufficient entries in eSum: " << hEnergy->GetEntries() << " for group " << ig << " energy point " << ie);
489 histFile->Close();
490 B2INFO("closed histFile; quitting");
491 return c_NotEnoughData;
492 }
493
494 //..Rebin if the maximum bin has too few entries
495 const double minPeakValue = 300.;
496 while (hEnergy->GetMaximum() < minPeakValue) {hEnergy->Rebin(2);}
497
498 //------------------------------------------------------------------------------
499 //..Find 50% range for Novo fit to find peak value. Rebin if there are too many bins
500 double fitFraction = 0.5;
501 int iLo = 1;
502 int iHi = hEnergy->GetNbinsX();
503 const int maxBinsToFit = 59;
504 while ((iHi - iLo) > maxBinsToFit) {
505 std::vector<int> iBins = eclNOptimalFitRange(hEnergy, fitFraction);
506 iLo = iBins[0];
507 iHi = iBins[1];
508 if ((iHi - iLo) > maxBinsToFit) {hEnergy->Rebin(2);}
509 }
510
511 //..Fit range
512 const int nBinsX = hEnergy->GetNbinsX();
513 const double fitlow = hEnergy->GetBinLowEdge(iLo);
514 const double fithigh = hEnergy->GetBinLowEdge(iHi + 1);
515 double sum = hEnergy->Integral(iLo, iHi);
516 B2INFO("group " << ig << " ie " << ie << " name " << name << " newName " << newName << " nBinsX " << nBinsX << " 50% target: " <<
517 fitFraction * hEnergy->GetEntries() << " iLo: " << iLo << " iHi: " << iHi << " sum: " << sum);
518
519 //------------------------------------------------------------------------------
520 //..Set up the fit. Initial values
521 double normalization = hEnergy->GetMaximum();
522 const int iMax = hEnergy->GetMaximumBin();
523 double peak = hEnergy->GetBinLowEdge(iMax);
524 double effSigma = 0.5 * (fithigh - fitlow);
525 double eta = 0.4; // typical value for good fits
526
527 func->SetParameters(normalization, peak, effSigma, eta);
528 func->SetParNames("normalization", "peak", "effSigma", "eta");
529
530 //..Parameter ranges
531 func->SetParLimits(1, fitlow, fithigh);
532 func->SetParLimits(2, 0., 2.*effSigma);
533 func->SetParLimits(3, -1, 3);
534
535 //..If there is only one crystal, fix eta to a special value
536 if (nCrysSumToFit == 1) {
537 const double etaForOneCrystal = 0.95;
538 func->SetParameter(3, etaForOneCrystal);
539 func->SetParLimits(3, etaForOneCrystal, etaForOneCrystal);
540 }
541
542 //..Fit
543 hEnergy->Fit(func, "LIEQ", "", fitlow, fithigh);
544
545 //------------------------------------------------------------------------------
546 //..Find the bias = sum of CalDigits minus sum of MC truth for this nCrys
547 // Bias is the mid point of the minimum range that contains 68.3% of events
548 double bias = 0.;
549 double biasSigma = 0.;
550 if (nCrysSumToFit < nCrysBins) {
551 TH1D* hBias = (TH1D*)biasSum->ProjectionY("hBias", nCrysSumToFit, nCrysSumToFit);
552 fitFraction = 0.683;
553 std::vector<int> jBins = eclNOptimalFitRange(hBias, fitFraction);
554 const double lowEdge = hBias->GetBinLowEdge(jBins[0]);
555 const double highEdge = hBias->GetBinLowEdge(jBins[1] + 1);
556 bias = 0.5 * (lowEdge + highEdge);
557 biasSigma = 0.5 * (highEdge - lowEdge);
558 }
559
560 //------------------------------------------------------------------------------
561 //..Resolution is half the smallest range that contains 68.3% of events.
562 // Speed up the process by only looking among the range that contains 75%.
563 fitFraction = 0.75;
564 std::vector<int> iBins = eclNOptimalFitRange(hEnergy, fitFraction);
565 const double resolution = eclNOptimalResolution(hEnergy, iBins[0], iBins[1]);
566
567 //..Correct for the bias from background when calculating fractional resolution.
568 // Note that variable resolution = (E_hi - E_lo)/E_gen is dimensionless,
569 // variable peak = E_peak/E_gen is dimensionless, but variable bias has units GeV
570 peak = func->GetParameter(1);
571 const double fractionalResolution = resolution / (peak - bias / genEnergy);
572
573 //..Record information if this is the best resolution (and is not the existing
574 // resolution case).
575 if (fractionalResolution < bestFractionalResolution and nCrysSumToFit != nCrysBins) {
576 bestFractionalResolution = fractionalResolution;
577 nOpt = nCrysSumToFit;
578 eSumOpt = hEnergy;
579 nFitBins = 1 + iHi - iLo;
580 biasForNOpt = bias;
581 biasSigmaForNOpt = biasSigma;
582 } else if (nCrysSumToFit == nCrysBins) {
583 existingFractionalResolution = fractionalResolution;
584 }
585
586 //------------------------------------------------------------------------------
587 //..Logic to decide what to do next
588
589 //..keep checking different N until resolution is this much worse than best value
590 const double resTolerance = 1.05;
591 if (nCrysSumToFit == initialnCrysSumToFit) {
592
593 //..After testing the previous nOptimal, try one fewer if possible
594 if (nCrysSumToFit > 1) {
595 nCrysSumToFit--;
596 } else {
597 nCrysSumToFit++;
598 }
599 } else if (nCrysSumToFit < initialnCrysSumToFit) {
600
601 //..Trying fewer crystals. If this is the best so far, try one
602 // fewer, if possible. Otherwise, try more crystals.
603 if (fractionalResolution < resTolerance* bestFractionalResolution and nCrysSumToFit > 1) {
604 nCrysSumToFit--;
605 } else if (initialnCrysSumToFit != nCrysMax) {
606 nCrysSumToFit = initialnCrysSumToFit + 1;
607 } else {
608 nCrysSumToFit = nCrysBins;
609 }
610 } else if (nCrysSumToFit < nCrysBins) { // nCrysSumToFit is always > initialnCrysSumToFit
611
612 //..Trying more crystals. If this is the best, try one more, if
613 // possible. Otherwise, do the current reconstruction case (nCrysSumToFit = nCrysBins)
614 if (fractionalResolution < resTolerance * bestFractionalResolution and nCrysSumToFit < nCrysMax) {
615 nCrysSumToFit++;
616 } else {
617 nCrysSumToFit = nCrysBins;
618 }
619 } else if (nCrysSumToFit == nCrysBins) {
620
621 //..This is the current resolution case, so we are done
622 nCrysSumToFit = 0;
623 }
624
625 } // while(nCrysSumToFit > 0)
626 B2INFO(" ig: " << ig << " ie: " << ie << " initial nOpt: " << initialnCrysSumToFit << " final nOpt: " << nOpt);
627
628 //-----------------------------------------------------------------------------------
629 //..Store everything in diagnostic histograms
630
631 //..Extract the function from the nOptimal histogram
632 TF1* funcOpt = (TF1*)eSumOpt->GetFunction("eclNOptimalNovo");
633
634 nOptimalPerGroup->SetBinContent(ig + 1, ie + 1, nOpt);
635
636 changeInNOptimal->SetBinContent(ig + 1, ie + 1, nOpt - initialnCrysSumToFit);
637 changeInNOptimal->SetBinError(ig + 1, ie + 1, 0.);
638
639 biasPerGroup->SetBinContent(ig + 1, ie + 1, biasForNOpt);
640 biasPerGroup->SetBinError(ig + 1, ie + 1, 0.);
641
642 sigmaBiasPerGroup->SetBinContent(ig + 1, ie + 1, biasSigmaForNOpt);
643 sigmaBiasPerGroup->SetBinError(ig + 1, ie + 1, 0.);
644
645 binsInFit->SetBinContent(ig + 1, ie + 1, nFitBins);
646 binsInFit->SetBinError(ig + 1, ie + 1, 0.);
647
648 maxEntriesPerHist->SetBinContent(ig + 1, ie + 1, eSumOpt->GetMaximum());
649 maxEntriesPerHist->SetBinError(ig + 1, ie + 1, 0.);
650
651 const double peakOpt = funcOpt->GetParameter(1);
652 const double peakOptUnc = funcOpt->GetParError(1);
653 peakPerGroup->SetBinContent(ig + 1, ie + 1, peakOpt);
654 peakPerGroup->SetBinError(ig + 1, ie + 1, peakOptUnc);
655
656 const double genE = generatedEPerGroup->GetBinContent(ig + 1, ie + 1);
657 const double peakCor = peakOpt - biasForNOpt / genE;
658 fracContainedEPerGroup->SetBinContent(ig + 1, ie + 1, peakCor);
659 fracContainedEPerGroup->SetBinError(ig + 1, ie + 1, peakOptUnc);
660
661 effSigmaPerGroup->SetBinContent(ig + 1, ie + 1, funcOpt->GetParameter(2));
662 effSigmaPerGroup->SetBinError(ig + 1, ie + 1, funcOpt->GetParError(2));
663
664 etaPerGroup->SetBinContent(ig + 1, ie + 1, funcOpt->GetParameter(3));
665 etaPerGroup->SetBinError(ig + 1, ie + 1, funcOpt->GetParError(3));
666
667 fitProbPerGroup->SetBinContent(ig + 1, ie + 1, funcOpt->GetProb());
668 fitProbPerGroup->SetBinError(ig + 1, ie + 1, 0.);
669
670 resolutionPerGroup->SetBinContent(ig + 1, ie + 1, bestFractionalResolution);
671 resolutionPerGroup->SetBinError(ig + 1, ie + 1, 0.);
672
673 existingResolutionPerGroup->SetBinContent(ig + 1, ie + 1, existingFractionalResolution);
674 existingResolutionPerGroup->SetBinError(ig + 1, ie + 1, 0.);
675
676 //..Write out to disk
677 histFile->cd();
678 eSumOpt->Write();
679
680 }
681 }
682
683 //-----------------------------------------------------------------------------------
684 //..Write out fit summary histograms
685 histFile->cd();
686 changeInNOptimal->Write();
687 biasPerGroup->Write();
688 sigmaBiasPerGroup->Write();
689 binsInFit->Write();
690 maxEntriesPerHist->Write();
691 peakPerGroup->Write();
692 effSigmaPerGroup->Write();
693 etaPerGroup->Write();
694 fitProbPerGroup->Write();
695 resolutionPerGroup->Write();
696 existingResolutionPerGroup->Write();
697 fracContainedEPerGroup->Write();
698
699 B2INFO("Wrote fit summary histograms");
700
701 //-----------------------------------------------------------------------------------
702 //-----------------------------------------------------------------------------------
703 //..Find bias and fraction contained energy in nOpt crystals for the energy bins
704 // adjacent to each of the group / energy points.
705
706 //..Histograms to store results
707 TH2F* fracContainedAdjacent[2];
708 TH2F* biasAdjacent[2];
709 const TString updown[2] = {"lower", "higher"};
710 for (int ia = 0; ia < 2; ia++) {
711 TString name = "fracContainedAdjacent_" + std::to_string(ia);
712 TString title = "peak fraction of energy contained in nOpt crystals, " + updown[ia] + " E;group;energy point";
713 fracContainedAdjacent[ia] = new TH2F(name, title, nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
714
715 name = "biasAdjacent_" + std::to_string(ia);
716 title = "bias (GeV) in nOpt crystals, " + updown[ia] + " E;group;energy point";
717 biasAdjacent[ia] = new TH2F(name, title, nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
718 }
719
720 //-----------------------------------------------------------------------------------
721 //..Loop over all groups and energy points
722 for (int ig = 0; ig < nCrystalGroups; ig++) {
723 for (int ie = 0; ie < nEnergies; ie++) {
724
725 //..nOpt for this point
726 const int nOpt = nOptimalPerGroup->GetBinContent(ig + 1, ie + 1);
727
728 //-----------------------------------------------------------------------------------
729 //..Two adjacent energy points
730 for (int ia = 0; ia < 2; ia++) {
731 const int ieAdj = ie + 2 * ia - 1; // ie +/- 1
732 double eFracPeakAdj = fracContainedEPerGroup->GetBinContent(ig + 1, ie + 1);
733 double biasAdj = biasPerGroup->GetBinContent(ig + 1, ie + 1);
734
735 //..Need this to be a valid energy point to do anything more
736 if (ieAdj >= 0 and ieAdj < nEnergies) {
737
738 //----------------------------------------------------------------------------
739 //..Find the bias
740 std::string biasName = "biasSum_" + std::to_string(ig) + "_" + std::to_string(ieAdj);
741 auto biasSum = getObjectPtr<TH2F>(biasName);
742 TH1D* hBias = (TH1D*)biasSum->ProjectionY("hBias", nOpt, nOpt);
743
744 //..Bias is the mid-point of the range containing 68.3% of events
745 double fitFraction = 0.683;
746 std::vector<int> jBins = eclNOptimalFitRange(hBias, fitFraction);
747 const double lowEdge = hBias->GetBinLowEdge(jBins[0]);
748 const double highEdge = hBias->GetBinLowEdge(jBins[1] + 1);
749 biasAdj = 0.5 * (lowEdge + highEdge);
750
751 //----------------------------------------------------------------------------
752 //..Get the eSum distribution to be fit
753 std::string name = "eSum_" + std::to_string(ig) + "_" + std::to_string(ieAdj);
754 auto eSum = getObjectPtr<TH2F>(name);
755 TH1D* hEnergy = (TH1D*)eSum->ProjectionY("hEnergy", nOpt, nOpt);
756 TString newName = name + "_" + std::to_string(nOpt);
757 hEnergy->SetName(newName);
758
759 B2INFO("fit adj ig = " << ig << " ie = " << ie << " ia = " << ia << " " << newName << " " << " entries = " << hEnergy->GetEntries()
760 << " GetMaxiumum = " << hEnergy->GetMaximum());
761
762 //..Rebin if the maximum bin has too few entries
763 const double minPeakValue = 300.;
764 while (hEnergy->GetMaximum() < minPeakValue) {hEnergy->Rebin(2);}
765
766 //..Find 50% range for Novo fit. Rebin if there are too many bins
767 fitFraction = 0.5;
768 int iLo = 1;
769 int iHi = hEnergy->GetNbinsX();
770 const int maxBinsToFit = 59;
771 while ((iHi - iLo) > maxBinsToFit) {
772 std::vector<int> iBins = eclNOptimalFitRange(hEnergy, fitFraction);
773 iLo = iBins[0];
774 iHi = iBins[1];
775 if ((iHi - iLo) > maxBinsToFit) {hEnergy->Rebin(2);}
776 }
777
778 //..Fit parameters
779 const double fitlow = hEnergy->GetBinLowEdge(iLo);
780 const double fithigh = hEnergy->GetBinLowEdge(iHi + 1);
781 double normalization = hEnergy->GetMaximum();
782 const int iMax = hEnergy->GetMaximumBin();
783 double peak = hEnergy->GetBinLowEdge(iMax);
784 double effSigma = 0.5 * (fithigh - fitlow);
785 double eta = 0.4; // typical value for good fits
786
787 func->SetParameters(normalization, peak, effSigma, eta);
788 func->SetParNames("normalization", "peak", "effSigma", "eta");
789
790 //..Parameter ranges
791 func->SetParLimits(1, fitlow, fithigh);
792 func->SetParLimits(2, 0., 2.*effSigma);
793 func->SetParLimits(3, -1, 3);
794
795 //..If there is only one crystal, fix eta to a special value
796 if (nOpt == 1) {
797 const double etaForOneCrystal = 0.95;
798 func->SetParameter(3, etaForOneCrystal);
799 func->SetParLimits(3, etaForOneCrystal, etaForOneCrystal);
800 }
801
802 //..Fit
803 hEnergy->Fit(func, "LIEQ", "", fitlow, fithigh);
804
805 //..Fraction of contained energy = peak of fit corrected for bias
806 const double peakAdj = func->GetParameter(1);
807 const double genE = generatedEPerGroup->GetBinContent(ig + 1, ieAdj + 1);
808 eFracPeakAdj = peakAdj - biasAdj / genE;
809 }
810
811 //..Store
812 fracContainedAdjacent[ia]->SetBinContent(ig + 1, ie + 1, eFracPeakAdj);
813 biasAdjacent[ia]->SetBinContent(ig + 1, ie + 1, biasAdj);
814 } // loop over adjacent energy points
815
816 B2INFO(" ig " << ig << " ie " << ie << " eFrac (3): "
817 << fracContainedAdjacent[0]->GetBinContent(ig + 1, ie + 1) << " "
818 << fracContainedEPerGroup->GetBinContent(ig + 1, ie + 1) << " "
819 << fracContainedAdjacent[1]->GetBinContent(ig + 1, ie + 1)
820 << " biases (3): "
821 << biasAdjacent[0]->GetBinContent(ig + 1, ie + 1) << " "
822 << biasPerGroup->GetBinContent(ig + 1, ie + 1) << " "
823 << biasAdjacent[1]->GetBinContent(ig + 1, ie + 1)
824 );
825 } // loop over ie
826 } // loop over ig
827
828 //-----------------------------------------------------------------------------------
829 //..Write out histograms with results for adjacent energy points
830 histFile->cd();
831 fracContainedAdjacent[0]->Write();
832 fracContainedAdjacent[1]->Write();
833 biasAdjacent[0]->Write();
834 biasAdjacent[1]->Write();
835
836 B2INFO("Wrote adjacent energy point summary histograms");
837
838
839 //-----------------------------------------------------------------------------------
840 //-----------------------------------------------------------------------------------
841 //..Prepare the payload contents
842 //..Upper energy boundaries are the mid-point of the log energies
843 auto boundaryE = create2Dvector<float>(nLeakReg, nEnergies);
844 for (int ireg = 0; ireg < nLeakReg; ireg++) {
845 B2INFO("Generated energies and boundaries for region = " << ireg);
846 for (int ie = 0; ie < nEnergies - 1; ie++) {
847 double logE0 = log(generatedE[ireg][ie]);
848 double logE1 = log(generatedE[ireg][ie + 1]);
849 double logBoundary = 0.5 * (logE0 + logE1);
850 boundaryE[ireg][ie] = exp(logBoundary);
851 B2INFO(" " << ie << " " << generatedE[ireg][ie] << " " << boundaryE[ireg][ie] << " GeV");
852 }
853
854 //..Last boundary is an arbitrarily large number
855 boundaryE[ireg][nEnergies - 1] = 9999.;
856 B2INFO(" " << nEnergies - 1 << " " << generatedE[ireg][nEnergies - 1] << " " << boundaryE[ireg][nEnergies - 1] << " GeV");
857 }
858
859 //..Group number of each cellID
860 std::vector<int> groupNumber;
861 for (int cellID = 1; cellID <= ECLElementNumbers::c_NCrystals; cellID++) {
862 const int iGroup = (int)(0.5 + groupNumberOfEachCellID->GetBinContent(cellID));
863 groupNumber.push_back(iGroup);
864 }
865
866 //..Copy results for each group into the payload histogram
867 TH2F nOptimal2D("nOptimal2D", "Optimal nCrys, energy bin vs group number", nCrystalGroups, 0, nCrystalGroups, nEnergies, 0.,
868 nEnergies);
869 for (int ig = 0; ig < nCrystalGroups; ig++) {
870 for (int ie = 0; ie < nEnergies; ie++) {
871 double nOpt = nOptimalPerGroup->GetBinContent(ig + 1, ie + 1);
872 nOptimal2D.SetBinContent(ig + 1, ie + 1, nOpt);
873 }
874 }
875
876 //..2D histograms for peak and bias corrections.
877 // Note that these have 3 x bins per group = nominal energy plus adjacent energies.
878 TH2F peakFracEnergy("peakFracEnergy", "peak contained energy over generated E, energy bin vs group number", 3 * nCrystalGroups, 0,
879 nCrystalGroups, nEnergies, 0., nEnergies);
880 TH2F bias("bias", "median bias (GeV), energy bin vs group number", 3 * nCrystalGroups, 0, nCrystalGroups, nEnergies, 0., nEnergies);
881 TH2F logPeakEnergy("logPeakEnergy", "log peak Energy (GeV), energy bin vs group number", 3 * nCrystalGroups, 0, nCrystalGroups,
882 nEnergies, 0., nEnergies);
883 for (int ig = 0; ig < nCrystalGroups; ig++) {
884 for (int ie = 0; ie < nEnergies; ie++) {
885 const int iy = ie + 1;
886 const int ix = 1 + 3 * ig; // three bins per group in payload histograms
887
888 //..Peak fractional energy and log peak contained energy (GeV)
889
890 //..Nominal ig and and ie
891 double peakAdj = fracContainedEPerGroup->GetBinContent(ig + 1, ie + 1);
892 const double genE = generatedEPerGroup->GetBinContent(ig + 1, ie + 1);
893 double logE = log(genE * peakAdj);
894
895 peakFracEnergy.SetBinContent(ix, iy, peakAdj);
896 logPeakEnergy.SetBinContent(ix, iy, logE);
897
898 //..Lower E adjacent point
899 double peakAdj0 = peakAdj;
900 double logE0 = logE;
901 if (ie > 1) {
902 peakAdj0 = fracContainedAdjacent[0]->GetBinContent(ig + 1, ie + 1);
903 const double genE0 = generatedEPerGroup->GetBinContent(ig + 1, ie);
904 logE0 = log(genE0 * peakAdj0);
905 }
906
907 peakFracEnergy.SetBinContent(ix + 1, iy, peakAdj0);
908 logPeakEnergy.SetBinContent(ix + 1, iy, logE0);
909
910 //..Higher E adjacent point
911 double peakAdj1 = peakAdj;
912 double logE1 = logE;
913 if (ie < nEnergies - 1) {
914 peakAdj1 = fracContainedAdjacent[1]->GetBinContent(ig + 1, ie + 1);
915 const double genE1 = generatedEPerGroup->GetBinContent(ig + 1, ie + 2);
916 logE1 = log(genE1 * peakAdj1);
917 }
918
919 peakFracEnergy.SetBinContent(ix + 2, iy, peakAdj1);
920 logPeakEnergy.SetBinContent(ix + 2, iy, logE1);
921
922 //..Bias
923 double medianBias = biasPerGroup->GetBinContent(ig + 1, ie + 1);
924 bias.SetBinContent(ix, iy, medianBias);
925 medianBias = biasAdjacent[0]->GetBinContent(ig + 1, ie + 1);
926 bias.SetBinContent(ix + 1, iy, medianBias);
927 medianBias = biasAdjacent[1]->GetBinContent(ig + 1, ie + 1);
928 bias.SetBinContent(ix + 2, iy, medianBias);
929 }
930 }
931
932 //..Store the histograms
933 histFile->cd();
934 nOptimal2D.Write();
935 peakFracEnergy.Write();
936 logPeakEnergy.Write();
937 bias.Write();
938 histFile->Close();
939
940 //-----------------------------------------------------------------------------------
941 //..The payload itself
942 ECLnOptimal* optimalNCrys = new ECLnOptimal();
943
944 //..Energy boundaries
945 for (int ireg = 0; ireg < nLeakReg; ireg++) {
946 std::vector<float> eUpperBoundaries;
947 for (int ie = 0; ie < nEnergies; ie++) {eUpperBoundaries.push_back(boundaryE[ireg][ie]);}
948 if (ireg == 0) {
949 optimalNCrys->setUpperBoundariesFwd(eUpperBoundaries);
950 } else if (ireg == 1) {
951 optimalNCrys->setUpperBoundariesBrl(eUpperBoundaries);
952 } else {
953 optimalNCrys->setUpperBoundariesBwd(eUpperBoundaries);
954 }
955 }
956
957 //..Group number of each cellID
958 optimalNCrys->setGroupNumber(groupNumber);
959
960 //..nOptimal histogram
961 optimalNCrys->setNOptimal(nOptimal2D);
962
963 //..Peak and bias correction histograms
964 optimalNCrys->setPeakFracEnergy(peakFracEnergy);
965 optimalNCrys->setBias(bias);
966 optimalNCrys->setLogPeakEnergy(logPeakEnergy);
967
968 //..Save the calibration
969 saveCalibration(optimalNCrys, "ECLnOptimal");
970 B2RESULT("eclNOptimalAlgorithm: successfully stored payload ECLnOptimal");
971
972 //..Done
973 return c_OK;
974}
975
976
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)
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
EResult
The result of calibration.
@ c_OK
Finished successfully =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Singleton class to cache database objects.
Definition: DBStore.h:31
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:54
void setInitializeActive(bool active)
Setter for m_initializeActive.
Definition: DataStore.cc:94
virtual EResult calibrate() override
..Run algorithm on events
DB object to store the optimal number of crystals to be used in a cluster energy sum,...
Definition: ECLnOptimal.h:23
void setNOptimal(const TH2F &nOptimal)
Set the 2D histogram of nOptimal.
Definition: ECLnOptimal.h:92
void setLogPeakEnergy(const TH2F &logPeakEnergy)
Set the 2D histogram of log of peak contained energy in GeV.
Definition: ECLnOptimal.h:101
void setGroupNumber(const std::vector< int > &groupNumber)
Set the vector of group numbers for each crystal.
Definition: ECLnOptimal.h:89
void setUpperBoundariesFwd(const std::vector< float > &eUpperBoundariesFwd)
Set vector of energies that specify energy ranges in the forward endcap.
Definition: ECLnOptimal.h:80
void setBias(const TH2F &bias)
Set the 2D histogram of beam background bias (reconstructed - true) (GeV)
Definition: ECLnOptimal.h:98
void setUpperBoundariesBrl(const std::vector< float > &eUpperBoundariesBrl)
Set vector of energies that specify energy ranges in the barrel.
Definition: ECLnOptimal.h:83
void setUpperBoundariesBwd(const std::vector< float > &eUpperBoundariesBwd)
Set vector of energies that specify energy ranges in the backward endcap.
Definition: ECLnOptimal.h:86
void setPeakFracEnergy(const TH2F &peakFracEnergy)
Set the 2D histogram of peak fractional contained energy.
Definition: ECLnOptimal.h:95
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
bool construct(Args &&... params)
Construct an object of type T in this StoreObjPtr, using the provided constructor arguments.
Definition: StoreObjPtr.h:119
static DBStore & Instance()
Instance of a singleton DBStore.
Definition: DBStore.cc:28
void updateEvent()
Updates all intra-run dependent objects.
Definition: DBStore.cc:142
void update()
Updates all objects that are outside their interval of validity.
Definition: DBStore.cc:79
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
bool isForward(int cellId)
Check whether the crystal is in forward ECL.
const int c_NCrystals
Number of crystals.
bool isBackward(int cellId)
Check whether the crystal is in backward ECL.
Abstract base class for different kinds of events.
STL namespace.