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