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