Belle II Software development
CDCDedxCosLayerAlgorithm.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/CDCDedxCosLayerAlgorithm.h>
10
11#include <TF1.h>
12#include <TLine.h>
13#include <TCanvas.h>
14#include <TH1I.h>
15#include <vector>
16
17using namespace Belle2;
18//-----------------------------------------------------------------
19// Implementation
20//-----------------------------------------------------------------
22 CalibrationAlgorithm("CDCDedxElectronCollector"),
23 isMethodSep(true),
24 isMakePlots(true),
25 isMerge(true),
26 isFixTrunc(false),
27 isUseTrunc(false),
28 m_truncMin(0.05),
29 m_truncMax(0.75),
30 m_cosBin(100),
31 m_cosMin(-1.0),
32 m_cosMax(1.0),
33 m_dedxBin(250),
34 m_dedxMin(0.0),
35 m_dedxMax(5.0),
36 m_suffix("")
37{
38 // Set module properties
39 setDescription("A calibration algorithm for CDC dE/dx electron cos(theta) dependence");
40
41}
42
43//-----------------------------------------------------------------
44// Run the calibration
45//-----------------------------------------------------------------
47{
48
50
51 if (!m_DBCosineCor.isValid())
52 B2FATAL("There is no valid previous payload for CDCDedxCosineCor");
53
54 B2INFO("Preparing dE/dx calibration for CDC dE/dx electron saturation");
55
56 // Get data objects
57 auto ttree = getObjectPtr<TTree>("tree");
58 if (ttree->GetEntries() < 100)return c_NotEnoughData;
59
60 std::vector<double>* lDedx = nullptr;
61 std::vector<int>* lLayer = nullptr;
62 double costh = 0.0;
63 int charge = 0;
64
65 ttree->SetBranchAddress("ldedx", &lDedx);
66 ttree->SetBranchAddress("lLayer", &lLayer);
67 ttree->SetBranchAddress("costh", &costh);
68 ttree->SetBranchAddress("charge", &charge);
69
70 // fill histograms, bin size may be arbitrary
71 TH1D* hCosth_neg = defineCosthHist("neg", "e-");
72 TH1D* hCosth_pos = defineCosthHist("pos", "e+");
73 TH1D* hCosth_all = defineCosthHist("all", "e-,e+");
74
75 std::array<std::vector<double>, m_kNGroups> corrFactor;
76
77 for (int ig = 0; ig < m_kNGroups; ig++) {
78 corrFactor[ig].assign(m_cosBin, 1.0);
79 }
80
81 for (int iter = 0; iter < 3; iter++) {
82 // make histograms to store dE/dx values in bins of cos(theta)
83 // bin size can be arbitrary, but for now just make uniform bins
84 std::array<std::vector<TH1D*>, m_kNGroups> hDedxCos_neg;
85 std::array<std::vector<TH1D*>, m_kNGroups> hDedxCos_pos;
86 std::array<std::vector<TH1D*>, m_kNGroups> hDedxCos_all;
87
88 defineHisto(hDedxCos_neg, Form("neg_iter%d", iter), "e-");
89 defineHisto(hDedxCos_pos, Form("pos_iter%d", iter), "e+");
90 defineHisto(hDedxCos_all, Form("all_iter%d", iter), "e-,e+");
91
92 std::array<TH1D*, m_kNGroups> hDedxGroup{};
93 for (int il = 0; il < m_kNGroups; il++) {
94 std::string title = Form("dedxhit dist (%s); dedxhit;entries", m_label[il].data());
95 hDedxGroup[il] = new TH1D(Form("hDedxGroup_%s_%s_iter%d", m_label[il].data(), m_suffix.data(), iter), "", m_dedxBin, m_dedxMin,
96 m_dedxMax);
97 hDedxGroup[il]->SetTitle(title.c_str());
98 }
99
100 const double binW = (m_cosMax - m_cosMin) / m_cosBin;
101
102 for (int i = 0; i < ttree->GetEntries(); ++i) {
103
104 ttree->GetEvent(i);
105
106 if (!lDedx || !lLayer) continue;
107 if (lDedx->size() != lLayer->size()) continue;
108 if (charge == 0) continue;
109
110 if (costh < TMath::Cos(150 * TMath::DegToRad()) ||
111 costh > TMath::Cos(17 * TMath::DegToRad())) continue;
112
113 int bin = int((costh - m_cosMin) / binW);
114 if (bin < 0 || bin >= int(m_cosBin)) continue;
115
116 for (size_t j = 0; j < lDedx->size(); ++j) {
117
118 double val = lDedx->at(j);
119 int lay = lLayer->at(j);
120
121 if (val <= 0) continue;
122
123 int ig = (lay < 8) ? 0 : ((lay < 14) ? 1 : 2);
124
125 val /= corrFactor[ig][bin];
126
127 if (isMethodSep) {
128 if (charge < 0)
129 hDedxCos_neg[ig][bin]->Fill(val);
130 else if (charge > 0)
131 hDedxCos_pos[ig][bin]->Fill(val);
132
133 }
134 hDedxCos_all[ig][bin]->Fill(val);
135
136 hDedxGroup[ig]->Fill(val);
137
138 }
139 if (iter == 0) {
140 // costh histo
141 if (isMethodSep) {
142 if (charge < 0) hCosth_neg->Fill(costh);
143 else if (charge > 0) hCosth_pos->Fill(costh);
144 }
145 hCosth_all->Fill(costh);
146 }
147 }
148
149 std::array<std::vector<double>, m_kNGroups> cosine;
150
151 for (int il = 0; il < m_kNGroups; ++il) {
152
153 int minGroup = 0, maxGroup = 0;
154 std::vector<double> vmean_neg, vmean_pos;
155
156 if (isFixTrunc && isUseTrunc) {
157 getTruncatedBins(hDedxGroup[il], minGroup, maxGroup);
158 hDedxGroup[il]->SetTitle(
159 Form("%s;%d;%d", hDedxGroup[il]->GetTitle(), minGroup, maxGroup));
160 }
161
162 cosine[il].reserve(m_cosBin);
163
164 for (unsigned int ibin = 0; ibin < m_cosBin; ++ibin) {
165
166 double mean = 1.0;
167
168 if (isMethodSep) {
169
170 double mean_neg = extractCosMean(hDedxCos_neg[il][ibin], minGroup, maxGroup);
171 double mean_pos = extractCosMean(hDedxCos_pos[il][ibin], minGroup, maxGroup);
172
173 bool has_neg = (hDedxCos_neg[il][ibin]->Integral() > 0);
174 bool has_pos = (hDedxCos_pos[il][ibin]->Integral() > 0);
175
176 if (has_neg && has_pos) mean = 0.5 * (mean_neg + mean_pos);
177 else if (has_neg) mean = mean_neg;
178 else if (has_pos) mean = mean_pos;
179 vmean_neg.push_back(mean_neg);
180 vmean_pos.push_back(mean_pos);
181
182 } else {
183 mean = extractCosMean(hDedxCos_all[il][ibin], minGroup, maxGroup);
184 }
185
186 cosine[il].push_back(mean);
187 if (mean > 0) corrFactor[il][ibin] *= (mean / 1.25);
188 }
189 if (isMethodSep) {
190 std::array<std::vector<double>, 3> cosMeanSets = {vmean_neg, vmean_pos, cosine[il]};
191 plotmeanChargeOverlay(cosMeanSets, m_label[il], iter);
192 }
193 }
194
195
196 if (isMakePlots) {
197
198 //1. dE/dx dist. for cosine bins
199
200 if (isMethodSep) {
201 plotdedxHist(hDedxCos_neg, Form("neg_iter%d", iter));
202 plotdedxHist(hDedxCos_pos, Form("pos_iter%d", iter));
203 }
204 plotdedxHist(hDedxCos_all, Form("all_iter%d", iter));
205
206 //3. Inner and Outer layer dE/dx distributions
207 plotLayerDist(hDedxGroup, iter);
208
209 //5. draw the relative constants
210 plotRelConst(corrFactor, iter);
211 }
212
213 }
214
215 for (int il = 0; il < m_kNGroups; ++il)
216 m_coscors.push_back(corrFactor[il]);
217
219
220 if (isMakePlots) {
221 //4. costh distribution
222 plotQaPars(hCosth_all, hCosth_pos, hCosth_neg);
223
224 //6. draw the final constants
226
227 //7. plot statistics related plots here
229
230 }
231
232 m_suffix.clear();
233
234 return c_OK;
235}
236
237//--------------------------------------------------
239{
240
241 int cruns = 0;
242 for (auto expRun : getRunList()) {
243 if (cruns == 0) B2INFO("CDCDedxBadWires: start exp " << expRun.first << " and run " << expRun.second << "");
244 cruns++;
245 }
246
247 const auto erStart = getRunList()[0];
248 int estart = erStart.first;
249 int rstart = erStart.second;
250
251 const auto erEnd = getRunList()[cruns - 1];
252 int eend = erEnd.first;
253 int rend = erEnd.second;
254
255 updateDBObjPtrs(1, rstart, estart);
256
257 m_runExp = Form("Range (%d:%d,%d:%d)", estart, rstart, eend, rend);
258 if (m_suffix.length() > 0) m_suffix = Form("%s_e%d_r%dr%d", m_suffix.data(), estart, rstart, rend);
259 else m_suffix = Form("e%d_r%dr%d", estart, rstart, rend);
260}
261
262//--------------------------------------------------
263TH1D* CDCDedxCosLayerAlgorithm::defineCosthHist(const std::string& tag, const std::string& chargeLabel)
264{
265
266 TH1D* hist = new TH1D(Form("hCosth_%s_%s", tag.c_str(), m_suffix.data()), " ", m_cosBin, m_cosMin, m_cosMax);
267 hist->SetTitle(Form("cos(#theta) dist (%s); cos(#theta); Entries", chargeLabel.c_str()));
268
269 return hist;
270}
271
272//--------------------------------------------------
273void CDCDedxCosLayerAlgorithm::defineHisto(std::array<std::vector<TH1D*>, m_kNGroups>& hdedx, const std::string& tag,
274 const std::string& chargeLabel)
275{
276
277 const double binW = (m_cosMax - m_cosMin) / m_cosBin;
278
279 for (int il = 0; il < m_kNGroups; il++) {
280 hdedx[il].reserve(m_cosBin);
281
282 for (unsigned int i = 0; i < m_cosBin; ++i) {
283 double coslow = i * binW + m_cosMin;
284 double coshigh = coslow + binW;
285
286 hdedx[il].push_back(new TH1D(Form("hDedxCos_%s_g%d_bin%d_%s", tag.c_str(), il, i, m_suffix.data()),
288
289 hdedx[il][i]->SetTitle(Form("%s dE/dx dist (%s) in costh (%0.02f, %0.02f)",
290 m_label[il].c_str(), chargeLabel.c_str(), coslow, coshigh));
291
292 // hdedx[il][i]->GetXaxis()->SetTitle("dE/dx");
293 // hdedx[il][i]->GetYaxis()->SetTitle("Entries");
294 }
295 }
296}
297
298//--------------------------------------------------
300{
301
302 for (unsigned int il = 0; il < m_kNGroups; il++) {
303 if (isMerge) {
304 unsigned int nbins = m_DBCosineCor->getSize(getRepresentativeLayer(il));
305 if (nbins != m_cosBin)
306 B2ERROR("merging failed because of unmatch bins (old " << m_cosBin << " new " << nbins << ")");
307
308 for (unsigned int ibin = 0; ibin < nbins; ibin++) {
309 double prev = m_DBCosineCor->getMean(getRepresentativeLayer(il), ibin);
310 B2INFO("Cosine Corr for " << m_label[il] << " Bin # " << ibin << ", Previous = " << prev << ", Relative = " << m_coscors[il][ibin]
311 << ", Merged = " << prev * m_coscors[il][ibin]);
312 m_coscors[il][ibin] *= prev;
313
314 }
315 }
316 }
317
318 //Saving constants
319 B2INFO("dE/dx calibration done for CDC dE/dx _eltron saturation");
320
321 std::vector<unsigned int> layerToGroup(56);
322
323 for (unsigned int layer = 0; layer < 56; layer++) {
324 if (layer < 8) layerToGroup[layer] = 0; // SL0
325 else if (layer < 14) layerToGroup[layer] = 1; // SL1
326 else layerToGroup[layer] = 2; // SL2-8
327 }
328
329 CDCDedxCosineCor* gain = new CDCDedxCosineCor(m_coscors, layerToGroup);
330 saveCalibration(gain, "CDCDedxCosineCor");
331}
332
333//--------------------------------------------------
334void CDCDedxCosLayerAlgorithm::plotdedxHist(std::array<std::vector<TH1D*>, 3>& hdedx, const std::string& tag)
335{
336
337 TCanvas ctmp("tmp", "tmp", 1200, 1200);
338 int nx = 2;
339 int ny = 2;
340 unsigned int nPads = nx * ny;
341
342 ctmp.Divide(nx, ny);
343 std::stringstream psname;
344
345 psname << Form("cdcdedx_coscorr_ldedx_%s_%s.pdf[", tag.c_str(), m_suffix.data());
346 ctmp.Print(psname.str().c_str());
347 psname.str("");
348 psname << Form("cdcdedx_coscorr_ldedx_%s_%s.pdf", tag.c_str(), m_suffix.data());
349
350 for (int il = 0; il < m_kNGroups; il++) {
351
352 for (unsigned int ic = 0; ic < m_cosBin; ic++) {
353
354 ctmp.cd(ic % nPads + 1);
355 hdedx[il][ic]->SetFillColor(4 + il);
356
357 hdedx[il][ic]->SetTitle(Form("%s;ldedx;entries", hdedx[il][ic]->GetTitle()));
358 hdedx[il][ic]->DrawClone("hist");
359
360 if (ic % nPads == nPads - 1 || ic == m_cosBin - 1) {
361 ctmp.Print(psname.str().c_str());
362 gPad->Clear("D");
363 ctmp.Clear("D");
364 }
365 }
366 }
367 psname.str("");
368 psname << Form("cdcdedx_coscorr_ldedx_%s_%s.pdf]", tag.c_str(), m_suffix.data());
369 ctmp.Print(psname.str().c_str());
370}
371
372//--------------------------------------------------
373void CDCDedxCosLayerAlgorithm::plotLayerDist(std::array<TH1D*, 3>& hDedxGroup, int iter)
374{
375
376 TCanvas cdedxlayer(Form("layerdedxhit_iter%d", iter), "Inner and Outer Layer dedxhit dist", 2400, 800);
377 cdedxlayer.Divide(3, 1);
378
379 for (int il = 0; il < m_kNGroups; il++) {
380 int minlay = 0, maxlay = 0;
381 if (isFixTrunc) {
382 minlay = std::stoi(hDedxGroup[il]->GetXaxis()->GetTitle());
383 maxlay = std::stoi(hDedxGroup[il]->GetYaxis()->GetTitle());
384 double lowedge = hDedxGroup[il]->GetXaxis()->GetBinLowEdge(minlay);
385 double upedge = hDedxGroup[il]->GetXaxis()->GetBinUpEdge(maxlay);
386 hDedxGroup[il]->SetTitle(Form("%s, trunc #rightarrow: %0.02f - %0.02f;dedxhit;entries", hDedxGroup[il]->GetTitle(), lowedge,
387 upedge));
388 }
389
390 cdedxlayer.cd(il + 1);
391 hDedxGroup[il]->SetFillColor(kYellow);
392 hDedxGroup[il]->Draw("histo");
393
394 if (isFixTrunc) {
395 TH1D* hDedxGroupC = (TH1D*)hDedxGroup[il]->Clone(Form("hDedxGroupC%d", il));
396 hDedxGroupC->GetXaxis()->SetRange(minlay, maxlay);
397 hDedxGroupC->SetFillColor(kAzure + 1);
398 hDedxGroupC->Draw("same histo");
399 }
400 }
401
402 cdedxlayer.SaveAs(Form("cdcdedx_coscorr_dedxlay%s_iter%d.pdf", m_suffix.data(), iter));
403 cdedxlayer.SaveAs(Form("cdcdedx_coscorr_dedxlay%s_iter%d.root", m_suffix.data(), iter));
404}
405
406//--------------------------------------------------
407void CDCDedxCosLayerAlgorithm::plotQaPars(TH1D* hCosth_all, TH1D* hCosth_pos, TH1D* hCosth_neg)
408{
409
410 TCanvas ceadist("ceadist", "Cosine distributions", 800, 600);
411 ceadist.cd();
412
413 TLegend* leg = new TLegend(0.6, 0.7, 0.8, 0.9);
414
415 // Always draw ALL first
416 if (hCosth_all) {
417 hCosth_all->SetFillColor(kYellow);
418 hCosth_all->SetLineColor(kBlack);
419 hCosth_all->SetStats(0);
420 hCosth_all->Draw("hist");
421 leg->AddEntry(hCosth_all, "all", "f");
422 }
423
424 // If method separation, overlay pos/neg
425 if (isMethodSep) {
426
427 if (hCosth_pos) {
428 hCosth_pos->SetLineColor(kRed);
429 hCosth_pos->SetFillStyle(0);
430 hCosth_pos->SetStats(0);
431 hCosth_pos->Draw("hist same");
432 leg->AddEntry(hCosth_pos, "pos", "l");
433 }
434
435 if (hCosth_neg) {
436 hCosth_neg->SetLineColor(kBlue);
437 hCosth_neg->SetFillStyle(0);
438 hCosth_neg->SetStats(0);
439 hCosth_neg->Draw("hist same");
440 leg->AddEntry(hCosth_neg, "neg", "l");
441 }
442 }
443
444 leg->Draw();
445
446 ceadist.SaveAs(Form("cdcdedx_coscorr_cosine_%s.pdf", m_suffix.data()));
447 ceadist.SaveAs(Form("cdcdedx_coscorr_cosine_%s.root", m_suffix.data()));
448}
449
450//--------------------------------------------------
451void CDCDedxCosLayerAlgorithm::plotRelConst(const std::array<std::vector<double>, m_kNGroups>& cosine, int iter)
452{
453 TCanvas cconst("cconst", "calibration Constants", 800, 600);
454 cconst.cd();
455
456 TLegend* leg = new TLegend(0.6, 0.8, 0.9, 0.9);
457 leg->SetBorderSize(0);
458 leg->SetFillStyle(0);
459
460 std::vector<TH1D*> hists;
461 std::vector<int> colors = {kRed, kBlue, kGreen + 2};
462 double ymax = 0.0;
463
464 for (int il = 0; il < m_kNGroups; il++) {
465
466 TH1D* h = new TH1D(Form("hconst_%d_%s", il, m_suffix.data()), "Relative constants; cos(#theta); constant", m_cosBin, m_cosMin,
467 m_cosMax);
468
469 // fill histogram
470 for (unsigned int jea = 0; jea < m_cosBin; jea++) {
471 if (jea < cosine[il].size())
472 h->SetBinContent(jea + 1, cosine[il].at(jea));
473 }
474
475 double hmax = h->GetMaximum();
476 if (hmax > ymax) ymax = hmax;
477
478
479 hists.push_back(h);
480 }
481
482 for (int il = 0; il < m_kNGroups; il++) {
483
484 hists[il]->SetLineColor(colors[il]);
485 hists[il]->SetStats(0);
486 hists[il]->SetLineWidth(2);
487
488 // draw
489 if (il == 0) {
490 hists[il]->SetMaximum(ymax + 0.01);
491 hists[il]->Draw("hist");
492 } else {
493 hists[il]->Draw("hist same");
494 }
495
496 leg->AddEntry(hists[il], m_label[il].data(), "l");
497 }
498
499 leg->Draw();
500
501 cconst.SaveAs(Form("cdcdedx_coscorr_relconst_%s_iter%d.pdf", m_suffix.data(), iter));
502 cconst.SaveAs(Form("cdcdedx_coscorr_relconst_%s_iter%d.root", m_suffix.data(), iter));
503
504 // cleanup
505 for (auto h : hists) delete h;
506}
507
508//--------------------------------------------------
510{
511
512 for (int il = 0; il < m_kNGroups; il++) {
513
514 unsigned int nbins = m_DBCosineCor->getSize(getRepresentativeLayer(il));
515
516 // --- Create histograms ---
517 TH1D* hnew = new TH1D(Form("hnew_%s", m_label[il].data()), Form("Final const: %s", m_label[il].data()), m_cosBin, m_cosMin,
518 m_cosMax);
519
520 TH1D* hold = new TH1D(Form("hold_%s", m_label[il].data()), Form("Final const: %s", m_label[il].data()), m_cosBin, m_cosMin,
521 m_cosMax);
522
523 for (unsigned int iea = 0; iea < nbins; iea++) {
524 double oldv = m_DBCosineCor->getMean(getRepresentativeLayer(il), iea);
525 double newv = m_coscors[il][iea];
526
527 hold->SetBinContent(iea + 1, oldv);
528 hnew->SetBinContent(iea + 1, newv);
529 }
530
531 // --- Ratio ---
532 TH1D* hratio = (TH1D*)hnew->Clone(Form("hratio_%s", m_label[il].data()));
533 hratio->Divide(hold);
534
535 TCanvas c(Form("c_%s", m_label[il].data()), Form("Final constants %s", m_label[il].data()), 1000, 500);
536 c.Divide(2, 1);
537 c.cd(1);
538
539 hnew->SetLineColor(kBlack);
540 hnew->SetStats(0);
541 hold->SetLineColor(kRed);
542 hold->SetStats(0);
543
544 double min = std::min(hnew->GetMinimum(), hold->GetMinimum());
545 double max = std::max(hnew->GetMaximum(), hold->GetMaximum());
546 hnew->GetYaxis()->SetRangeUser(min * 0.95, max * 1.05);
547
548 hnew->Draw("hist");
549 hold->Draw("hist same");
550
551 auto leg = new TLegend(0.6, 0.75, 0.85, 0.88);
552 leg->SetBorderSize(0);
553 leg->SetFillStyle(0);
554 leg->AddEntry(hnew, "New", "l");
555 leg->AddEntry(hold, "Old", "l");
556 leg->Draw();
557
558 c.cd(2);
559
560 hratio->SetLineColor(kBlue);
561 hratio->SetStats(0);
562 hratio->SetTitle(Form("Ratio: new/old, %s;cos(#theta); New / Old", m_label[il].data()));
563 hratio->GetYaxis()->SetRangeUser(0.2, 1.2);
564 hratio->Draw("hist");
565
566 TLine* line = new TLine(m_cosMin, 1.0, m_cosMax, 1.0);
567 line->SetLineStyle(2);
568 line->Draw();
569
570 c.SaveAs(Form("cdcdedx_coscorr_fconsts_%s_%s.pdf", m_label[il].data(), m_suffix.data()));
571 c.SaveAs(Form("cdcdedx_coscorr_fconsts_%s_%s.root", m_label[il].data(), m_suffix.data()));
572
573 // cleanup
574 delete hnew;
575 delete hold;
576 delete hratio;
577 delete line;
578 }
579}
580
581//--------------------------------------------------
582void CDCDedxCosLayerAlgorithm::plotmeanChargeOverlay(const std::array<std::vector<double>, 3>& mean, const std::string& sltag,
583 int iter)
584{
585 TCanvas cconst(Form("cconst_%s_iter%d", sltag.c_str(), iter), "calibration Constants", 800, 600);
586 cconst.cd();
587
588 TLegend* leg = new TLegend(0.6, 0.8, 0.9, 0.9);
589 leg->SetBorderSize(0);
590 leg->SetFillStyle(0);
591
592 std::vector<TH1D*> hists;
593
594 std::array<std::string, 3> labels = {"e^{+}", "e^{-}", "Average"};
595 std::vector<int> colors = {kRed, kBlue, kBlack + 2};
596
597 for (int il = 0; il < 3; il++) {
598
599 TH1D* h = new TH1D(Form("hconst_%d_%s_%s_iter%d", il, sltag.c_str(), m_suffix.data(), iter),
600 Form("Relative mean %s, iter %d; cos(#theta); constant", sltag.c_str(), iter),
602
603 for (unsigned int jea = 0; jea < m_cosBin; jea++) {
604 if (jea < mean[il].size())
605 h->SetBinContent(jea + 1, mean[il].at(jea));
606 }
607
608 h->SetLineColor(colors[il]);
609 h->SetStats(0);
610 h->SetLineWidth(2);
611
612 if (il == 0) h->Draw("hist");
613 else h->Draw("hist same");
614
615 leg->AddEntry(h, labels[il].c_str(), "l");
616 hists.push_back(h);
617 }
618
619 leg->Draw();
620
621 cconst.SaveAs(Form("cdcdedx_coscorr_relmean_iter%d_%s_%s.pdf", iter, sltag.c_str(), m_suffix.data()));
622 cconst.SaveAs(Form("cdcdedx_coscorr_relmean_iter%d_%s_%s.root", iter, sltag.c_str(), m_suffix.data()));
623
624 for (auto h : hists) delete h;
625 delete leg;
626}
627
628//------------------------------------
630{
631
632 TCanvas cstats("cstats", "cstats", 1000, 500);
633 cstats.SetBatch(kTRUE);
634 cstats.Divide(2, 1);
635
636 cstats.cd(1);
637 auto hestats = getObjectPtr<TH1I>("hestats");
638 if (hestats) {
639 hestats->SetName(Form("hestats_%s", m_suffix.data()));
640 hestats->SetStats(0);
641 hestats->DrawCopy("");
642 }
643
644 cstats.cd(2);
645 auto htstats = getObjectPtr<TH1I>("htstats");
646 if (htstats) {
647 htstats->SetName(Form("htstats_%s", m_suffix.data()));
648 htstats->SetStats(0);
649 htstats->DrawCopy("");
650 }
651 cstats.Print(Form("cdcdedx_coscorr_stats_%s.pdf", m_suffix.data()));
652}
653
654//--------------------------------------------------
655void CDCDedxCosLayerAlgorithm::getTruncatedBins(TH1D* hist, int& binlow, int& binhigh)
656{
657
658 //calculating truncation average
659 double sum = hist->Integral();
660 if (sum <= 0 || hist->GetNbinsX() <= 0) {
661 binlow = 1; binhigh = 1;
662 return ;
663 }
664
665 binlow = 1.0; binhigh = 1.0;
666 double sumPer5 = 0.0, sumPer75 = 0.0;
667 for (int ibin = 1; ibin <= hist->GetNbinsX(); ibin++) {
668 double bcdedx = hist->GetBinContent(ibin);
669 if (sumPer5 <= m_truncMin * sum) {
670 sumPer5 += bcdedx;
671 binlow = ibin;
672 }
673 if (sumPer75 <= m_truncMax * sum) {
674 sumPer75 += bcdedx;
675 binhigh = ibin;
676 }
677 }
678 return;
679}
680
681//--------------------------
682double CDCDedxCosLayerAlgorithm::getTruncationMean(TH1D* hist, int binlow, int binhigh)
683{
684
685 //calculating truncation average
686 if (hist->Integral() < 100) return 1.0;
687
688 if (binlow <= 0 || binhigh > hist->GetNbinsX())return 1.0;
689
690 double binweights = 0., sumofbc = 0.;
691 for (int ibin = binlow; ibin <= binhigh; ibin++) {
692 double bcdedx = hist->GetBinContent(ibin);
693 if (bcdedx > 0) {
694 binweights += (bcdedx * hist->GetBinCenter(ibin));
695 sumofbc += bcdedx;
696 }
697 }
698 if (sumofbc > 0) return binweights / sumofbc;
699 else return 1.0;
700}
701
702double CDCDedxCosLayerAlgorithm::extractCosMean(TH1D*& hist, int fixedLow, int fixedHigh)
703{
704 if (!hist || hist->Integral() <= 0) return 1.0;
705
706 // Default for SL0, SL1: simple mean
707 if (!isUseTrunc) {
708
709 hist->SetTitle(Form("%s, mean = %0.5f", hist->GetTitle(), hist->GetMean()));
710 return hist->GetMean();
711 }
712
713 int minbin = 1, maxbin = 1;
714 if (isFixTrunc) {
715 minbin = fixedLow;
716 maxbin = fixedHigh;
717 } else {
718 getTruncatedBins(hist, minbin, maxbin);
719 }
720
721 double mean = getTruncationMean(hist, minbin, maxbin);
722
723 hist->SetTitle(Form("%s, mean = %0.5f;%d;%d", hist->GetTitle(), mean, minbin, maxbin));
724
725 return mean;
726}
void plotdedxHist(std::array< std::vector< TH1D * >, m_kNGroups > &hdedx, const std::string &tag)
function to draw the dE/dx histogram in costh bins
void plotRelConst(const std::array< std::vector< double >, m_kNGroups > &cosine, int iter)
Plot relative calibration constants vs costh for all SL groups (overlayed)
void plotQaPars(TH1D *hCosth_all, TH1D *hCosth_pos, TH1D *hCosth_neg)
function to costh distribution for Inner/Outer layer
double m_truncMax
upper threshold on truncation
void plotLayerDist(std::array< TH1D *, m_kNGroups > &hdedxlay, int iter)
function to draw dedx dist.
double m_truncMin
lower threshold on truncation
void getTruncatedBins(TH1D *hist, int &binlow, int &binhigh)
function to get bins of truncation from histogram
void defineHisto(std::array< std::vector< TH1D * >, m_kNGroups > &hdedx, const std::string &tag, const std::string &chargeLabel)
function to define dE/dx histograms
void getExpRunInfo()
function to get extract calibration run/exp
static constexpr int m_kNGroups
SL grouping: inner (SL0), middle (SL1), outer (SL2–8)
int m_dedxBin
number of bins for dedx histogram
double m_cosMax
max cosine angle for cal
void plotmeanChargeOverlay(const std::array< std::vector< double >, 3 > &cosine_pos, const std::string &sltag, int iter)
Plot overlay of positive, negative, and average cosine means for one SL group.
unsigned int getRepresentativeLayer(unsigned int igroup) const
Representative CDC layer for each SL group (used to access group-wise constants): SL0 => 1,...
CDCDedxCosLayerAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
TH1D * defineCosthHist(const std::string &tag, const std::string &chargeLabel)
function to define cosine histograms
std::string m_suffix
add suffix to all plot name
double extractCosMean(TH1D *&hist, int fixedLow=1, int fixedHigh=1)
Extract mean dE/dx vs costh for a given group from the histogram.
double getTruncationMean(TH1D *hist, int binlow, int binhigh)
function to get truncated mean
DBObjPtr< CDCDedxCosineCor > m_DBCosineCor
Electron saturation correction DB object.
void plotConstants()
function to draw the old/new final constants
bool isFixTrunc
true = fix window for all out/inner layers
std::string m_runExp
add run and exp to title of plot
const std::array< std::string, m_kNGroups > m_label
add inner/outer superlayer label
double m_cosMin
min cosine angle for cal
void plotEventStats()
function to draw the stats plots
virtual EResult calibrate() override
Cosine algorithm.
bool isMethodSep
if e+e- need to be consider sep
double m_dedxMax
max dedx range for gain cal
void createPayload()
function to store new payload after full calibration
bool isMakePlots
produce plots for status
bool isMerge
merge payload at the of calibration
bool isUseTrunc
true if truncated mean for SL0,1
double m_dedxMin
min dedx range for gain cal
unsigned int m_cosBin
number of bins across cosine range
std::vector< std::vector< double > > m_coscors
final vectors of calibration
dE/dx cosine gain calibration constants
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.