Belle II Software development
CDCDedxWireGainAlgorithm.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#include <reconstruction/calibration/CDCDedxWireGainAlgorithm.h>
10
11using namespace Belle2;
12using namespace CDC;
13using namespace std;
14
15//-----------------------------------------------------------------
16// Implementation
17//-----------------------------------------------------------------
19 CalibrationAlgorithm("CDCDedxElectronCollector"),
20 m_isMakePlots(true),
21 m_isMerge(true),
22 m_isWireTruc(false),
23 m_dedxBins(250),
24 m_dedxMin(0.0),
25 m_dedxMax(5.0),
26 m_truncMin(0.05),
27 m_truncMax(0.75),
28 m_suffix("")
29{
30 // Set module properties
31 setDescription("A calibration algorithm for CDC dE/dx wire gains");
32}
33
34//-----------------------------------------------------------------
35// Run the calibration
36//-----------------------------------------------------------------
38{
39
40 //if not checking then copy function here;
42
43 if (!m_DBBadWires.isValid() || !m_DBWireGains.isValid())
44 B2FATAL("There is no valid payload for BadWires and/or Wirgain");
45
46 // Get data objects
47 auto ttree = getObjectPtr<TTree>("tree");
48 if (ttree->GetEntries() < 100)return c_NotEnoughData;
49
50 vector<int>* wire = 0;
51 ttree->SetBranchAddress("wire", &wire);
52
53 vector<double>* dedxhit = 0;
54 ttree->SetBranchAddress("dedxhit", &dedxhit);
55
56 // dedxhit vector to store dE/dx values for each wire
57 map<int, vector<double>> wirededx;
58
59 //IL = Inner Layer and OL = Outer Layer
60 array<TH1D*, 2> hdedxL;
61 string label[2] = {"IL", "OL"};
62 for (int il = 0; il < 2; il++)
63 hdedxL[il] = new TH1D(Form("hdedx%s_%s", label[il].data(), m_suffix.data()), "", m_dedxBins, m_dedxMin, m_dedxMax);
64
65 for (int i = 0; i < ttree->GetEntries(); ++i) {
66 ttree->GetEvent(i);
67 for (unsigned int j = 0; j < wire->size(); ++j) {
68 int jwire = wire->at(j);
69 double jhitdedx = dedxhit->at(j);
70 wirededx[jwire].push_back(jhitdedx);
71 //wire # 1279 end of inner layers
72 if (jwire < 1280) hdedxL[0]->Fill(jhitdedx);
73 else hdedxL[1]->Fill(jhitdedx);
74 }
75 }
76
77 //return if ~10% low stats or dead wires
78 int minstat = 0;
79 for (unsigned int jw = 0; jw < c_nwireCDC; ++jw)
80 if (wirededx[jw].size() <= 100) minstat++;
81
82 if (minstat > 0.10 * c_nwireCDC) return c_NotEnoughData;
83
84 //25-75 is average bin # for truncation
85 array<unsigned int, 2> minbinL, maxbinL;
86 for (int il = 0; il < 2; il++) {
87 getTruncatedBins(hdedxL[il], minbinL[il], maxbinL[il]);
88 hdedxL[il]->SetTitle(Form("%s(%s);%d;%d", label[il].data(), m_suffix.data(), minbinL[il], maxbinL[il]));
89 }
90
91 vector<double> vrel_mean;
92 vector<double> vdedx_means;
93
94 vector<TH1D*> hdedxhit(c_nwireCDC);
95
96 B2INFO("Creating CDCGeometryPar object");
98
99 vector<double> layermean(c_maxNSenseLayers);
100
101 int activelayers = 0;
102 double layeravg = 0.0;
103
104 int jwire = -1;
105 for (unsigned int il = 0; il < c_maxNSenseLayers; ++il) {
106
107 int activewires = 0;
108 layermean[il] = 0.0;
109
110 for (unsigned int iw = 0; iw < cdcgeo.nWiresInLayer(il); ++iw) {
111
112 jwire++;
113 hdedxhit[jwire] = new TH1D(Form("h%s_w%d", m_suffix.data(), jwire), "", m_dedxBins, m_dedxMin, m_dedxMax);
114
115 for (unsigned int ih = 0; ih < wirededx[jwire].size(); ++ih) {
116 hdedxhit[jwire]->Fill(wirededx[jwire][ih]);
117 }
118
119 unsigned int minbin, maxbin;
120 if (!m_isWireTruc) {
121 getTruncatedBins(hdedxhit[jwire], minbin, maxbin);
122 } else {
123 if (jwire < 1280) {
124 minbin = minbinL[0];
125 maxbin = maxbinL[0];
126 } else {
127 minbin = minbinL[1];
128 maxbin = maxbinL[1];
129 }
130 }
131 hdedxhit[jwire]->SetTitle(Form("dedxhit-dist, wire: %d (%s);%d;%d", jwire, m_suffix.data(), minbin, maxbin));
132
133 double dedxmean;
134 if (m_DBBadWires->getBadWireStatus(jwire) == kTRUE) dedxmean = 0.0;
135
136 else dedxmean = getTruncationMean(hdedxhit[jwire], minbin, maxbin);
137 vrel_mean.push_back(dedxmean);
138
139 double prewg = m_DBWireGains->getWireGain(jwire);
140 if (prewg > 0.0 && m_isMerge) {
141 vdedx_means.push_back(dedxmean * prewg);
142 B2INFO("merged-wireGain: [" << jwire << "], prewgvious = " << prewg << ", rel = " << dedxmean << ", merged = " << vdedx_means.at(
143 jwire));
144 } else vdedx_means.push_back(dedxmean);
145
146 //calculate layer average for active wires
147 if (vdedx_means.at(jwire) > 0) {
148 layermean[il] += vdedx_means.at(jwire);
149 activewires++;
150 }
151 }
152
153 if (activewires > 0) layermean[il] /= activewires;
154 else layermean[il] = 1.0;
155
156 //calculate outer layer average for active layer
157 if (il >= 8 && layermean[il] > 0) {
158 layeravg += layermean[il];
159 activelayers++;
160 }
161
162 }
163
164 //normalisation of wiregains to get outergain ~1.0
165 if (activelayers > 0) layeravg /= activelayers;
166
167 for (unsigned int iw = 0; iw < c_nwireCDC; ++iw) {
168 vrel_mean.at(iw) /= layeravg;
169 vdedx_means.at(iw) /= layeravg;
170 }
171
172 if (m_isMakePlots) {
173
174 //1. Inner and Outer layer dE/dx distributions
175 plotLayerDist(hdedxL);
176
177 //2. wiredist -> 14336 (good, bad)
178 plotWireDist(hdedxhit, vrel_mean);
179
180 //3. wiregains draw (1D, dist)
181 plotWireGain(vdedx_means, vrel_mean, layeravg);
182
183 //4. layer gain plot
184 plotLayerGain(layermean, layeravg);
185
186 //5. wiregains per layer
187 plotWGPerLayer(vdedx_means, layermean, layeravg);
188
189 //6. Statstics plot
191 }
192
193 createPayload(vdedx_means);
194 m_suffix.clear();
195 return c_OK;
196}
197
198//------------------------------------
200{
201
202 int cruns = 0;
203 for (auto expRun : getRunList()) {
204 if (cruns == 0) B2INFO("CDCDedxWireGain: start exp " << expRun.first << " and run " << expRun.second << "");
205 cruns++;
206 }
207
208 const auto erStart = getRunList()[0];
209 int estart = erStart.first;
210 int rstart = erStart.second;
211
212 const auto erEnd = getRunList()[cruns - 1];
213 int rend = erEnd.second;
214
215 updateDBObjPtrs(1, rstart, estart);
216
217 if (m_suffix.length() > 0) m_suffix = Form("%s_e%d_r%dr%d", m_suffix.data(), estart, rstart, rend);
218 else m_suffix = Form("e%d_r%dr%d", estart, rstart, rend);
219}
220
221//--------------------------
222void CDCDedxWireGainAlgorithm::createPayload(const vector<double>& vdedx_means)
223{
224
225 B2INFO("dE/dx Calibration done for " << vdedx_means.size() << " CDC wires");
226 CDCDedxWireGain* gains = new CDCDedxWireGain(vdedx_means);
227 saveCalibration(gains, "CDCDedxWireGain");
228}
229
230//--------------------------
231double CDCDedxWireGainAlgorithm::getTruncationMean(TH1D* hdedxhit, int binlow, int binhigh)
232{
233
234 //calculating truncation average
235 if (hdedxhit->Integral() < 100) return 1.0;
236
237 if (binlow <= 0 || binhigh > hdedxhit->GetNbinsX()) return 1.0;
238
239 double binweights = 0., sumofbc = 0.;
240 for (int ibin = binlow; ibin <= binhigh; ibin++) {
241 double bcdedx = hdedxhit->GetBinContent(ibin);
242 if (bcdedx > 0) {
243 binweights += (bcdedx * hdedxhit->GetBinCenter(ibin));
244 sumofbc += bcdedx;
245 }
246 }
247 if (sumofbc > 0) return binweights / sumofbc;
248 else return 1.0;
249}
250
251//--------------------------
252void CDCDedxWireGainAlgorithm::getTruncatedBins(TH1D* hdedxhit, unsigned int& binlow, unsigned int& binhigh)
253{
254
255 //calculating truncation average
256 double sum = hdedxhit->Integral();
257 if (sum <= 0 || hdedxhit->GetNbinsX() <= 0) {
258 binlow = 1; binhigh = 1;
259 return ;
260 }
261
262 binlow = 1.0; binhigh = 1.0;
263 double sumPer5 = 0.0, sumPer75 = 0.0;
264 for (int ibin = 1; ibin <= hdedxhit->GetNbinsX(); ibin++) {
265 double bcdedx = hdedxhit->GetBinContent(ibin);
266 if (sumPer5 <= m_truncMin * sum) {
267 sumPer5 += bcdedx;
268 binlow = ibin;
269 }
270 if (sumPer75 <= m_truncMax * sum) {
271 sumPer75 += bcdedx;
272 binhigh = ibin;
273 }
274 }
275 return;
276}
277
278//--------------------------
279void CDCDedxWireGainAlgorithm::plotLayerDist(array<TH1D*, 2> hdedxL)
280{
281
282 TCanvas cldedx("cldedx", "IL/OL dedxhit dist", 900, 400);
283 cldedx.Divide(2, 1);
284
285 for (int il = 0; il < 2; il++) {
286 cldedx.cd(il + 1);
287
288 int minbin = stoi(hdedxL[il]->GetXaxis()->GetTitle());
289 int maxbin = stoi(hdedxL[il]->GetYaxis()->GetTitle());
290 double lowedge = hdedxL[il]->GetXaxis()->GetBinLowEdge(minbin);
291 double upedge = hdedxL[il]->GetXaxis()->GetBinUpEdge(maxbin);
292
293 hdedxL[il]->SetFillColor(kYellow);
294 hdedxL[il]->SetTitle(Form("%s, trunc(%0.02f - %0.02f);dedxhit;entries", hdedxL[il]->GetTitle(), lowedge, upedge));
295 hdedxL[il]->Draw("histo");
296
297 TH1D* hdedxLC = (TH1D*)hdedxL[il]->Clone(Form("%s_c", hdedxL[il]->GetName()));
298 hdedxLC->GetXaxis()->SetRange(minbin, maxbin);
299 hdedxLC->SetFillColor(kAzure + 1);
300 hdedxLC->Draw("same histo");
301 }
302
303 cldedx.SaveAs(Form("cdcdedx_wgcal_layerdedx_%s.pdf", m_suffix.data()));
304}
305
306//------------------------------------
307void CDCDedxWireGainAlgorithm::plotWireDist(const vector<TH1D*>& hist, const vector<double>& vrel_mean)
308{
309
310 TCanvas ctmp(Form("cdcdedx_%s", m_suffix.data()), "", 1200, 1200);
311 ctmp.Divide(4, 4);
312 ctmp.SetBatch(kTRUE);
313
314 stringstream psname;
315 psname << Form("cdcdedx_wgcal_%s.pdf[", m_suffix.data());
316 ctmp.Print(psname.str().c_str());
317 psname.str("");
318 psname << Form("cdcdedx_wgcal_%s.pdf", m_suffix.data());
319
320 for (unsigned int iw = 0; iw < hist.size(); iw++) {
321
322 int minbin = stoi(hist[iw]->GetXaxis()->GetTitle());
323 int maxbin = stoi(hist[iw]->GetYaxis()->GetTitle());
324
325 hist[iw]->SetFillColor(kYellow - 9);
326 hist[iw]->SetTitle(Form("%s, rel. #mu_{trunc} %0.03f;dedxhit;entries", hist[iw]->GetTitle(), vrel_mean.at(iw)));
327
328 if (m_DBBadWires->getBadWireStatus(iw) == kTRUE) {
329 hist[iw]->SetLineColor(kRed);
330 hist[iw]->SetLineWidth(2);
331 }
332 ctmp.cd(iw % 16 + 1);
333 hist[iw]->DrawCopy();
334
335 TH1D* hdedxhitC = (TH1D*)hist[iw]->Clone(Form("%sC", hist[iw]->GetName()));
336 hdedxhitC->GetXaxis()->SetRange(minbin, maxbin);
337 hdedxhitC->SetFillColor(kAzure + 1);
338 hdedxhitC->DrawCopy("same histo");
339
340 if (((iw + 1) % 16 == 0) || iw == (hist.size() - 1)) {
341 ctmp.Print(psname.str().c_str());
342 ctmp.Clear("D");
343 }
344 delete hist[iw];
345 delete hdedxhitC;
346 }
347
348 psname.str("");
349 psname << Form("cdcdedx_wgcal_%s.pdf]", m_suffix.data());
350 ctmp.Print(psname.str().c_str());
351}
352
353//------------------------------------
354void CDCDedxWireGainAlgorithm::plotWireGain(const vector<double>& vdedx_means, const vector<double>& vrel_mean, double layeravg)
355{
356
357 //saving final constants in a histograms for validation
358 TCanvas cwconst("cwconst", "", 900, 500);
359 TCanvas cwconstvar("cwconstvar", "", 500, 400);
360
361 array<TH1D*, 2> hconstpw, hconstpwvar;
362
363 for (int i = 0; i < 2; i++) {
364
365 hconstpw[i] = new TH1D(Form("hconstpw_%d_%s", i, m_suffix.data()), "", c_nwireCDC, -0.5, 14335.5);
366 hconstpw[i]->SetTitle(Form("wiregain const (%s); wire numbers;<dedxhit>", m_suffix.data()));
367 if (m_isMerge
368 && i == 0) hconstpw[i]->SetTitle(Form("merged wiregain rel-const (%s), avg = %0.03f; wire numbers;<dedxhit>", m_suffix.data(),
369 layeravg));
370
371 hconstpwvar[i] = new TH1D(Form("hconstpwvar_%s", m_suffix.data()), "", 400, -0.5, 2.5);
372 hconstpwvar[i]->SetTitle(Form("wiregain const (%s); wire gains; nentries", m_suffix.data()));
373 if (m_isMerge
374 && i == 0) hconstpwvar[i]->SetTitle(Form("merged wiregain rel-const (%s), avg = %0.03f; wire gains; nentries", m_suffix.data(),
375 layeravg));
376
377 for (unsigned int iw = 0; iw < c_nwireCDC; iw++) {
378
379 double gain = vdedx_means.at(iw);
380 if (m_isMerge && i == 1) gain = vrel_mean.at(iw);
381 hconstpw[i]->SetBinContent(iw + 1, gain);
382 hconstpwvar[i]->Fill(gain);
383
384 if (iw % 500 == 0) hconstpw[i]->GetXaxis()->SetBinLabel(iw + 1, Form("w%d", iw + 1));
385 }
386
387 hconstpw[i]->SetLineColor(i * 2 + 2);
388 hconstpw[i]->LabelsOption("u", "X");
389 hconstpw[i]->GetYaxis()->SetRangeUser(-0.1, 3.5);
390 hconstpw[i]->LabelsDeflate();
391
392 hconstpwvar[i]->SetFillColor(i * 2 + 2);
393 hconstpwvar[i]->Scale(1 / hconstpwvar[i]->GetMaximum());
394
395 }
396 cwconst.cd();
397 cwconst.SetGridy(1);
398 hconstpw[0]->Draw("");
399 if (m_isMerge) hconstpw[1]->Draw("same");
400
401 cwconstvar.cd();
402 hconstpwvar[0]->Draw("hist");
403 if (m_isMerge) hconstpwvar[1]->Draw("hist same");
404
405 cwconst.SaveAs(Form("cdcdedx_wgcal_wireconst_%s.pdf", m_suffix.data()));
406 cwconstvar.SaveAs(Form("cdcdedx_wgcal_wireconstvar_%s.pdf", m_suffix.data()));
407}
408
409//------------------------------------
410void CDCDedxWireGainAlgorithm::plotLayerGain(const vector<double>& layermean, double layeravg)
411{
412
413 TH1D hlayeravg(Form("hlayeravg_%s", m_suffix.data()), "", layermean.size(), -0.5, 55.5);
414 hlayeravg.SetTitle(Form("layer gain avg (%s); layer numbers;<dedxhit>", m_suffix.data()));
415
416 for (unsigned int il = 0; il < layermean.size(); il++) {
417 hlayeravg.SetBinContent(il + 1, layermean[il]);
418 if (il % 2 == 0 || il == layermean.size() - 1) hlayeravg.GetXaxis()->SetBinLabel(il + 1, Form("L%d", il));
419 }
420
421 TCanvas clayeravg("clayeravg", "clayeravg", 800, 500);
422 clayeravg.SetGridy(1);
423 clayeravg.cd();
424 gStyle->SetOptStat("ne");
425 hlayeravg.LabelsOption("u", "X");
426 hlayeravg.SetLineColor(kBlue);
427 hlayeravg.GetYaxis()->SetRangeUser(-0.1, 3.5);
428 hlayeravg.SetTitle(Form("%s, avg = %0.03f (abs)", hlayeravg.GetTitle(), layeravg));
429 hlayeravg.LabelsDeflate();
430 hlayeravg.Draw("");
431 TLine* tl = new TLine(-0.5, layeravg, 55.5, layeravg);
432 tl->SetLineColor(kRed);
433 tl->DrawClone("same");
434 clayeravg.SaveAs(Form("cdcdedx_wgcal_layeravg_%s.pdf", m_suffix.data()));
435 delete tl;
436}
437
438//------------------------------------
439void CDCDedxWireGainAlgorithm::plotWGPerLayer(const vector<double>& vdedx_means, const vector<double>& layermean, double layeravg)
440{
441
443
444 TCanvas clconst("clconst", "", 800, 500);
445 clconst.Divide(2, 2);
446 clconst.SetBatch(kTRUE);
447
448 stringstream psnameL;
449 psnameL << Form("cdcdedx_wgcal_layerconst_%s.pdf[", m_suffix.data());
450 clconst.Print(psnameL.str().c_str());
451 psnameL.str(""); psnameL << Form("cdcdedx_wgcal_layerconst_%s.pdf", m_suffix.data());
452
453 int jwire = 0;
454
455 for (unsigned int il = 0; il < layermean.size(); ++il) {
456
457 unsigned int nwires = cdcgeo.nWiresInLayer(il);
458 TH1D hconstpl(Form("hconstpwvar_l%d_%s", il, m_suffix.data()), "", nwires, jwire, jwire + nwires);
459 hconstpl.SetTitle(Form("abs-const, layer: %d (%s); wire numbers;<dedxhit>", il, m_suffix.data()));
460
461 for (unsigned int iw = 0; iw < nwires; ++iw) {
462 hconstpl.SetBinContent(iw + 1, vdedx_means.at(jwire));
463 if (il < 32) {
464 if (iw % 10 == 0) hconstpl.GetXaxis()->SetBinLabel(iw + 1, Form("w%d", jwire));
465 } else {
466 if (iw % 15 == 0) hconstpl.GetXaxis()->SetBinLabel(iw + 1, Form("w%d", jwire));
467 }
468 jwire++;
469 }
470
471 double lmean = layermean.at(il) / layeravg;
472
473 clconst.cd(il % 4 + 1);
474 gStyle->SetOptStat("ne");
475
476 hconstpl.SetTitle(Form("%s, avg = %0.03f", hconstpl.GetTitle(), lmean));
477
478 if (il < 8) hconstpl.GetYaxis()->SetRangeUser(-0.1, 4.0);
479 else hconstpl.GetYaxis()->SetRangeUser(-0.1, 2.0);
480 hconstpl.SetFillColor(kAzure - 1);
481 hconstpl.LabelsOption("u", "X");
482 hconstpl.DrawCopy("hist");
483
484 TLine* tlc = new TLine(jwire - nwires, lmean, jwire, lmean);
485 tlc->SetLineColor(kRed);
486 tlc->DrawClone("same");
487
488 if ((il + 1) % 4 == 0 || (il + 1) == layermean.size()) {
489 clconst.Print(psnameL.str().c_str());
490 clconst.Clear("D");
491 }
492
493 delete tlc;
494 }
495
496 psnameL.str("");
497 psnameL << Form("cdcdedx_wgcal_layerconst_%s.pdf]", m_suffix.data());
498 clconst.Print(psnameL.str().c_str());
499}
500
501//--------------------------
503{
504
505 TCanvas cstats("cstats", "cstats", 800, 400);
506 cstats.SetBatch(kTRUE);
507 cstats.Divide(2, 1);
508
509 cstats.cd(1);
510 auto hestats = getObjectPtr<TH1I>("hestats");
511 if (hestats) {
512 hestats->SetName(Form("hestats_%s", m_suffix.data()));
513 hestats->SetStats(0);
514 hestats->DrawCopy("");
515 }
516
517 cstats.cd(2);
518 auto htstats = getObjectPtr<TH1I>("htstats");
519 if (htstats) {
520 htstats->SetName(Form("htstats_%s", m_suffix.data()));
521 htstats->SetStats(0);
522 htstats->DrawCopy("");
523 }
524
525 cstats.Print(Form("cdcdedx_wgcal_stats_%s.pdf", m_suffix.data()));
526}
void plotWGPerLayer(const std::vector< double > &vdedx_means, const std::vector< double > &layermean, double layeravg)
function to draw WG per layer
void plotLayerGain(const std::vector< double > &layermean, double layeravg)
function to draw layer gains
void plotLayerDist(std::array< TH1D *, 2 > hdedxL)
function to draw dE/dx for inner/outer layer
bool m_isMerge
merge payload at the time of calibration
double m_truncMax
max trunc range for mean
double m_truncMin
min trunc range for mean
void getTruncatedBins(TH1D *hdedxhit, unsigned int &binlow, unsigned int &binhigh)
function to get bins of truncation from histogram
void plotWireGain(const std::vector< double > &vdedx_means, const std::vector< double > &vrel_mean, double layeravg)
function to draw wire gains
void getExpRunInfo()
function to get extract calibration run/exp
void createPayload(const std::vector< double > &vdedx_tmeans)
function to finally store new payload after full calibration
void plotWireDist(const std::vector< TH1D * > &hist, const std::vector< double > &vrel_mean)
function to draw dE/dx histograms for each wire
bool m_isMakePlots
Save arithmetic and truncated mean for the 'dedx' values.
DBObjPtr< CDCDedxBadWires > m_DBBadWires
Bad wire DB object.
std::string m_suffix
suffix string to separate plots
double getTruncationMean(TH1D *hdedxhit, int binlow, int binhigh)
function to get mean of truncation from histogram
int m_dedxBins
number of bins for dedx histogram
DBObjPtr< CDCGeometry > m_cdcGeo
Geometry of CDC.
CDCDedxWireGainAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
DBObjPtr< CDCDedxWireGain > m_DBWireGains
Wire gain DB object.
void plotEventStats()
function to draw statstics
virtual EResult calibrate() override
Wire gain algorithm.
double m_dedxMax
max dedx range for wiregain cal
double m_dedxMin
min dedx range for wiregain cal
bool m_isWireTruc
method of trunc range for mean
dE/dx wire gain calibration constants
The Class for CDC Geometry Parameters.
unsigned nWiresInLayer(int layerId) const
Returns wire numbers in a layer.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void updateDBObjPtrs(const unsigned int event, const int run, const int experiment)
Updates any DBObjPtrs by calling update(event) for DBStore.
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 successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
Abstract base class for different kinds of events.
STL namespace.