Belle II Software development
CDCDedx1DCellAlgorithm.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/CDCDedx1DCellAlgorithm.h>
10
11#include <TCanvas.h>
12#include <TLegend.h>
13#include <TMath.h>
14#include <TPad.h>
15#include <TRandom.h>
16#include <TStyle.h>
17
18#include <cmath>
19
20using namespace Belle2;
21
22//-----------------------------------------------------------------
23// Implementation
24//-----------------------------------------------------------------
26 CalibrationAlgorithm("CDCDedxElectronCollector"),
27 m_eaMin(-TMath::Pi() / 2),
28 m_eaMax(+TMath::Pi() / 2),
29 m_eaB(316),
30 m_dedxMin(0.0),
31 m_dedxMax(5.0),
32 m_dedxBin(250),
33 m_ptMax(8.0),
34 m_cosMax(1.0),
35 m_truncMin(0.05),
36 m_truncMax(0.75),
37 m_binSplit(3),
38 m_chargeType(0),
39 m_adjustFac(1.00),
40 isFixTrunc(false),
41 isVarBins(true),
42 isRotSymm(false),
43 isMakePlots(true),
44 isPrintLog(false),
45 isMerge(true),
46 m_suffix("")
47{
48 // Set module properties
49 setDescription("A calibration algorithm for the CDC dE/dx entrance angle cleanup correction");
50}
51
52//-----------------------------------------------------------------
53// Run the calibration
54//-----------------------------------------------------------------
56{
57
59
60 if (!m_DBOneDCell.isValid())
61 B2FATAL("There is no valid previous payload for CDCDedx1DCell");
62
63 //reading radiative electron collector TREE
64 auto ttree = getObjectPtr<TTree>("tree");
65 if (!ttree) return c_NotEnoughData;
66
67 std::vector<double>* dedxhit = 0, *enta = 0;
68 std::vector<int>* layer = 0;
69 double pt = 0, costh = 0;
70
71 ttree->SetBranchAddress("dedxhit", &dedxhit);
72 ttree->SetBranchAddress("entaRS", &enta);
73 ttree->SetBranchAddress("layer", &layer);
74 ttree->SetBranchAddress("pt", &pt);
75 ttree->SetBranchAddress("costh", &costh);
76
77 //repair nbins if they are not divisible accordingly
80
81 //Settings of variables bins
83
84 if (isPrintLog) {
85 B2INFO("Superlayer0 bins: " << m_eaBinLocal[0]);
86 B2INFO("Superlayer1 bins: " << m_eaBinLocal[1]);
87 B2INFO("Superlayers2-8 bins: " << m_eaBinLocal[2]);
88 }
89
90 // dedxhit vector to store dE/dx values for each enta bin
91 std::array<std::vector<TH1D*>, m_kNGroups> hdedxhit;
92 std::array<TH1D*, m_kNGroups> hdedxlay{};
93 std::array<TH1D*, m_kNGroups> hentalay{};
94
95 TH2D* hptcosth = new TH2D("hptcosth", "pt vs costh dist;pt;costh", 1000, -8.0, 8.0, 1000, -1.0, 1.0);
96
97 defineHisto(hdedxhit, hdedxlay, hentalay);
98
99 //Star filling histogram defined above
100 for (int i = 0; i < ttree->GetEntries(); ++i) {
101
102 ttree->GetEvent(i);
103
104 if (std::abs(costh) > m_cosMax) continue;
105
106 // remove wide angle bhabha tracks
107 // double mom = pt/sqrt(1-costh*costh);
108 // if(abs(pt)<2.4 && abs(mom)>3.6)continue;
109
110 if (std::abs(pt) > m_ptMax) continue;
111
112 //change to random 10%
113 int rand = gRandom->Integer(100);
114 if (rand < 10) hptcosth->Fill(pt, costh);
115
116 for (unsigned int j = 0; j < dedxhit->size(); ++j) {
117
118 if (dedxhit->at(j) == 0) continue;
119
120 double entaval = enta->at(j);
121 //Mapped bin corresponds to entaval
122 int ibin = std::floor((entaval - m_eaMin) / m_eaBW);
123 if (ibin < 0 || ibin >= m_eaBin) continue;
124
125 int lay = layer->at(j);
126 int mL = (lay < 8) ? 0 : ((lay < 14) ? 1 : 2);
127
128 hdedxlay[mL]->Fill(dedxhit->at(j));
129 if (rand < 10) hentalay[mL]->Fill(entaval);
130
131 int jbinea = ibin;
132 if (isVarBins) jbinea = m_binIndex[mL].at(ibin);
133 hdedxhit[mL][jbinea]->Fill(dedxhit->at(j));
134 }
135 }
136
137 for (int il = 0; il < m_kNGroups; il++) {
138
139 int minlay = 0, maxlay = 0;
140
141 if (isFixTrunc) {
142 getTruncatedBins(hdedxlay[il], minlay, maxlay);
143 hdedxlay[il]->SetTitle(Form("%s;%d;%d", hdedxlay[il]->GetTitle(), minlay, maxlay));
144 }
145
146 std::vector<double>tempconst;
147 tempconst.reserve(m_eaBinLocal[il]);
148
149 for (int iea = 0; iea < m_eaBinLocal[il]; iea++) {
150
151 int jea = iea;
152
153 // rotation symmtery for 1<->3 and 4<->2 but only symmetric bin
154 if (!isVarBins && isRotSymm) jea = rotationalBin(m_eaBinLocal[il], jea);
155
156 TH1D* htemp = (TH1D*)hdedxhit[il][jea]->Clone(Form("h_%s_b%d_c", m_label[il].data(), jea));
157
158 int minbin = 1, maxbin = 1;
159 if (isFixTrunc) {
160 minbin = minlay;
161 maxbin = maxlay;
162 } else {
163 //extract truncation window per bin
164 getTruncatedBins(htemp, minbin, maxbin);
165 }
166
167 double dedxmean;
168 dedxmean = getTruncationMean(htemp, minbin, maxbin);
169 tempconst.push_back(dedxmean);
170
171 hdedxhit[il][iea]->SetTitle(Form("%s, #mu_{truc} = %0.5f;%d;%d", hdedxhit[il][iea]->GetTitle(), dedxmean, minbin, maxbin));
172 }
173
174 //Expending constants
175 std::vector<double>layerconst;
176 layerconst.reserve(m_eaBin);
177
178 for (int iea = 0; iea < m_eaBin; iea++) {
179 int jea = iea;
180 if (isVarBins) jea = m_binIndex[il].at(iea);
181 layerconst.push_back(tempconst.at(jea));
182 }
183
184 // plot the rel constants var/sym bins
185 if (isMakePlots) plotRelConst(tempconst, layerconst, il);
186 m_onedcors.push_back(layerconst);
187
188 layerconst.clear();
189 tempconst.clear();
190 }
191
192 //Saving final constants
194
195 if (isMakePlots) {
196
197 //1. dE/dx dist. for entrance angle bins
198 plotdedxHist(hdedxhit);
199
200 //3. Inner and Outer layer dE/dx distributions
201 plotLayerDist(hdedxlay);
202
203 //4. entrance angle distribution sym/var bins
204 plotQaPars(hentalay, hptcosth);
205
206 //6. draw the final constants
208
209 //7. plot statistics related plots here
211 }
212
213 for (int il = 0; il < m_kNGroups; il++) {
214 delete hentalay[il];
215 delete hdedxlay[il];
216 for (int iea = 0; iea < m_eaBinLocal[il]; iea++)
217 delete hdedxhit[il][iea];
218 }
219
220 delete hptcosth;
221 m_eaBinLocal.clear();
222 for (int il = 0; il < m_kNGroups; il++) {
223 m_binValue[il].clear();
224 m_binIndex[il].clear();
225 }
226 m_suffix.clear();
227
228 return c_OK;
229}
230
231//--------------------------------------------------
233{
234
235 int cruns = 0;
236 for (auto expRun : getRunList()) {
237 if (cruns == 0) B2INFO("CDCDedxBadWires: start exp " << expRun.first << " and run " << expRun.second << "");
238 cruns++;
239 }
240
241 const auto erStart = getRunList()[0];
242 int estart = erStart.first;
243 int rstart = erStart.second;
244
245 const auto erEnd = getRunList()[cruns - 1];
246 int eend = erEnd.first;
247 int rend = erEnd.second;
248
249 updateDBObjPtrs(1, rstart, estart);
250
251 m_runExp = Form("Range (%d:%d,%d:%d)", estart, rstart, eend, rend);
252 if (m_suffix.length() > 0) m_suffix = Form("%s_e%d_r%dr%d", m_suffix.data(), estart, rstart, rend);
253 else m_suffix = Form("e%d_r%dr%d", estart, rstart, rend);
254}
255
256//--------------------------------------------------
258{
259
260 std::map<int, std::vector<double>> bounds;
261 std::map<int, std::vector<int>> steps;
262
263 // number of intervals for each unique configuration
264 const std::array<int, 2> nDev{8, 4};
265
266 // config 0 -> used for SL = 0,1
267 bounds[0] = {0, 108, 123, 133, 158, 183, 193, 208, 316};
268 steps[0] = {9, 3, 2, 1, 1, 2, 3, 9};
269
270 // config 1 -> used for SL = 2-8
271 bounds[1] = {0, 38, 158, 278, 316};
272 steps[1] = {2, 1, 1, 2};
273
274 // map SL -> configuration
275 const std::array<int, m_kNGroups> configIndex = {0, 0, 1};
276
277 for (int il = 0; il < m_kNGroups; il++) {
278
279 int icfg = configIndex[il];
280
281 std::vector<double> scaledBounds = bounds[icfg];
282 for (int ibin = 0; ibin <= nDev[icfg]; ibin++) {
283 scaledBounds[ibin] = scaledBounds[ibin] * m_binSplit;
284 }
285
286 int ieaprime = -1, temp = -99, ibin = 0;
287
288 double pastbin = m_eaMin;
289 m_binValue[il].push_back(pastbin);
290
291 for (int iea = 0; iea < m_eaBin; iea++) {
292
293 if (isVarBins) {
294 if (iea % int(scaledBounds[ibin + 1]) == 0 && iea > 0) ibin++;
295 int diff = iea - int(scaledBounds[ibin]);
296 if (diff % steps[icfg][ibin] == 0) ieaprime++;
297 } else {
298 ieaprime = iea;
299 }
300
301 m_binIndex[il].push_back(ieaprime);
302
303 if (ieaprime != temp) {
304 double binwidth = m_eaBW;
305 if (isVarBins) binwidth = m_eaBW * steps[icfg][ibin];
306 double binvalue = pastbin + binwidth;
307 pastbin = binvalue;
308 if (std::abs(binvalue) < 1e-5) binvalue = 0;
309 m_binValue[il].push_back(binvalue);
310 }
311 temp = ieaprime;
312 }
313
314 m_eaBinLocal.push_back(int(m_binValue[il].size()) - 1);
315 }
316
317 if (isMakePlots) plotMergeFactor(bounds, nDev, steps);
318}
319
320//--------------------------------------------------
321
322void CDCDedx1DCellAlgorithm::defineHisto(std::array<std::vector<TH1D*>, m_kNGroups>& hdedxhit,
323 std::array<TH1D*, m_kNGroups>& hdedxlay,
324 std::array<TH1D*, m_kNGroups>& hentalay)
325{
326 for (int il = 0; il < m_kNGroups; il++) {
327
328 std::string title = Form("dedxhit dist (%s): %s ; dedxhit;entries", m_label[il].data(), m_runExp.data());
329 hdedxlay[il] = new TH1D(Form("hdedxlay%s", m_label[il].data()), "", m_dedxBin, m_dedxMin, m_dedxMax);
330 hdedxlay[il]->SetTitle(title.c_str());
331
332 Double_t* nvarBins = m_binValue[il].data();
333
334 if (isVarBins)
335 title = Form("entaRS dist (variable bins): %s: (%s); entaRS (#alpha);entries", m_label[il].data(), m_runExp.data());
336 else
337 title = Form("entaRS dist (sym. bins): %s: (%s); entaRS (#alpha);entries", m_label[il].data(), m_runExp.data());
338
339 hentalay[il] = new TH1D(Form("hentalay%s", m_label[il].data()), "", m_eaBinLocal[il], nvarBins);
340 hentalay[il]->SetTitle(title.c_str());
341
342 for (int iea = 0; iea < m_eaBinLocal[il]; iea++) {
343
344 double min = m_binValue[il].at(iea);
345 double max = m_binValue[il].at(iea + 1);
346 double width = max - min;
347
348 if (isPrintLog)
349 B2INFO("bin: " << iea << " ], min:" << min << " , max: " << max << " , width: " << width);
350
351 title = Form("%s: entaRS = (%0.03f to %0.03f)", m_label[il].data(), min, max);
352
353 hdedxhit[il].push_back(new TH1D(Form("hdedxhit_%s_bin%d", m_label[il].data(), iea),
355
356 hdedxhit[il][iea]->SetTitle(title.c_str());
357 }
358 }
359}
360
361//--------------------------------------------------
362void CDCDedx1DCellAlgorithm::getTruncatedBins(TH1D* hist, int& binlow, int& binhigh)
363{
364
365 //calculating truncation average
366 double sum = hist->Integral();
367 if (sum <= 0 || hist->GetNbinsX() <= 0) {
368 binlow = 1; binhigh = 1;
369 return ;
370 }
371
372 binlow = 1.0; binhigh = 1.0;
373 double sumPer5 = 0.0, sumPer75 = 0.0;
374 for (int ibin = 1; ibin <= hist->GetNbinsX(); ibin++) {
375 double bcdedx = hist->GetBinContent(ibin);
376 if (sumPer5 <= m_truncMin * sum) {
377 sumPer5 += bcdedx;
378 binlow = ibin;
379 }
380 if (sumPer75 <= m_truncMax * sum) {
381 sumPer75 += bcdedx;
382 binhigh = ibin;
383 }
384 }
385 return;
386}
387
388//--------------------------
389double CDCDedx1DCellAlgorithm::getTruncationMean(TH1D* hist, int binlow, int binhigh)
390{
391
392 //calculating truncation average
393 if (hist->Integral() < 100) return 1.0;
394
395 if (binlow <= 0 || binhigh > hist->GetNbinsX())return 1.0;
396
397 double binweights = 0., sumofbc = 0.;
398 for (int ibin = binlow; ibin <= binhigh; ibin++) {
399 double bcdedx = hist->GetBinContent(ibin);
400 if (bcdedx > 0) {
401 binweights += (bcdedx * hist->GetBinCenter(ibin));
402 sumofbc += bcdedx;
403 }
404 }
405 if (sumofbc > 0) return binweights / sumofbc;
406 else return 1.0;
407}
408
409//--------------------------------------------------
411{
412
413 B2INFO("dE/dx one cell calibration: Generating payloads");
414
415 for (unsigned int il = 0; il < m_kNGroups; il++) {
416
417 if (isMerge) {
418 unsigned int nbins = m_DBOneDCell->getNBins(getRepresentativeLayer(il));
419
420 if (int(nbins) != m_eaBin)
421 B2ERROR("merging failed because of unmatch bins (old " << m_eaBin << " new " << nbins << ")");
422
423 for (unsigned int iea = 0; iea < nbins; iea++) {
424 double prev = m_DBOneDCell->getMean(getRepresentativeLayer(il), iea);
425 m_onedcors[il][iea] *= prev;
426 }
427 }
428
429 double binsize = TMath::Pi() / m_onedcors[il].size();
430
431 auto computeAverages = [&](double angLow, double angHigh) {
432 unsigned int binLow = std::floor((angLow + TMath::Pi() / 2.0) / binsize);
433 unsigned int binHigh = std::floor((angHigh + TMath::Pi() / 2.0) / binsize);
434 double sum_new = 0.0, sum_prev = 0.0;
435 int count = 0;
436
437 for (unsigned int iea = binLow; iea < binHigh; ++iea) {
438 sum_new += m_onedcors[il][iea];
439 sum_prev += m_DBOneDCell->getMean(getRepresentativeLayer(il), iea);
440 ++count;
441 }
442
443 double avg_new = (count > 0) ? sum_new / count : 1.0;
444 double avg_prev = (count > 0) ? sum_prev / count : 1.0;
445 return std::make_pair(avg_new, avg_prev);
446 };
447
448 double negLow = -0.5, negHigh = -0.2;
449 double posLow = 0.2, posHigh = 0.5;
450
451 if (il == 2) {
452 negLow = -0.75; negHigh = -0.25;
453 posLow = 0.25; posHigh = 0.75;
454 }
455
456 auto [avgNewNeg, avgPrevNeg] = computeAverages(negLow, negHigh);
457 auto [avgNewPos, avgPrevPos] = computeAverages(posLow, posHigh);
458
459 double avgNew = (avgNewNeg + avgNewPos) / 2.0;
460 double avgPrev = (avgPrevNeg + avgPrevPos) / 2.0;
461 double scaleFactor = avgPrev / avgNew;
462
463 for (int iea = 0; iea < m_eaBin; iea++) {
464 m_onedcors[il][iea] *= scaleFactor;
465 }
466
467 if (m_chargeType > 0)
468 for (int ie = 0; ie < m_eaBin / 2; ie++) m_onedcors[il][ie] *= m_adjustFac;
469 if (m_chargeType < 0)
470 for (int ie = m_eaBin / 2; ie < m_eaBin; ie++) m_onedcors[il][ie] *= m_adjustFac;
471
472 }
473
474 //Saving constants
475 B2INFO("dE/dx Calibration done for CDCDedx1DCell");
476 std::vector<unsigned int> layerToGroup(56);
477
478 for (unsigned int layer = 0; layer < 56; layer++) {
479 if (layer < 8) layerToGroup[layer] = 0; // SL0
480 else if (layer < 14) layerToGroup[layer] = 1; // SL1
481 else layerToGroup[layer] = 2; // SL2-8
482 }
483
484 CDCDedx1DCell* gain = new CDCDedx1DCell(0, m_onedcors, layerToGroup);
485 saveCalibration(gain, "CDCDedx1DCell");
486}
487
488//--------------------------------------------------
489void CDCDedx1DCellAlgorithm::plotMergeFactor(std::map<int, std::vector<double>> bounds,
490 const std::array<int, 2> nDev,
491 std::map<int, std::vector<int>> steps)
492{
493 TCanvas cmfactor("cmfactor", "Merging factors", 800, 400);
494 cmfactor.Divide(2, 1);
495
496 std::array<TH1I*, 2> hists{};
497
498 for (int icfg = 0; icfg < 2; icfg++) {
499 Double_t* nvarBins = bounds[icfg].data();
500
501 hists[icfg] = new TH1I(Form("hist_cfg%d", icfg), "", nDev[icfg], nvarBins);
502
503 if (icfg == 0)
504 hists[icfg]->SetTitle("Merging factor for SL0/SL1 bins;binindex;merge-factors");
505 else
506 hists[icfg]->SetTitle("Merging factor for SL2-8 bins;binindex;merge-factors");
507
508 for (int ibin = 0; ibin < nDev[icfg]; ibin++) {
509 hists[icfg]->SetBinContent(ibin + 1, steps[icfg][ibin]);
510 }
511
512 cmfactor.cd(icfg + 1);
513 hists[icfg]->SetFillColor(kYellow);
514 hists[icfg]->Draw("hist");
515 }
516
517 cmfactor.SaveAs(Form("cdcdedx_1dcell_mergefactor%s.pdf", m_suffix.data()));
518 cmfactor.SaveAs(Form("cdcdedx_1dcell_mergefactor%s.root", m_suffix.data()));
519
520 for (int icfg = 0; icfg < 2; icfg++) {
521 delete hists[icfg];
522 }
523}
524
525//--------------------------------------------------
526void CDCDedx1DCellAlgorithm::plotdedxHist(std::array<std::vector<TH1D*>, 3>& hdedxhit)
527{
528
529 TCanvas ctmp("tmp", "tmp", 1200, 1200);
530 ctmp.Divide(4, 4);
531 std::stringstream psname;
532
533 psname << Form("cdcdedx_1dcell_dedxhit%s.pdf[", m_suffix.data());
534 ctmp.Print(psname.str().c_str());
535 psname.str("");
536 psname << Form("cdcdedx_1dcell_dedxhit%s.pdf", m_suffix.data());
537
538 for (int il = 0; il < m_kNGroups; il++) {
539
540 for (int jea = 0; jea < m_eaBinLocal[il]; jea++) {
541
542 int minbin = std::stoi(hdedxhit[il][jea]->GetXaxis()->GetTitle());
543 int maxbin = std::stoi(hdedxhit[il][jea]->GetYaxis()->GetTitle());
544
545 ctmp.cd(jea % 16 + 1);
546 hdedxhit[il][jea]->SetFillColor(4 + il);
547
548 hdedxhit[il][jea]->SetTitle(Form("%s;dedxhit;entries", hdedxhit[il][jea]->GetTitle()));
549 hdedxhit[il][jea]->DrawClone("hist");
550 TH1D* htempC = (TH1D*)hdedxhit[il][jea]->Clone(Form("%sc2", hdedxhit[il][jea]->GetName()));
551 htempC->GetXaxis()->SetRange(minbin, maxbin);
552 htempC->SetFillColor(kGray);
553 htempC->DrawClone("same hist");
554
555 if (jea % 16 == 15 || (jea == m_eaBinLocal[il] - 1)) {
556 ctmp.Print(psname.str().c_str());
557 gPad->Clear("D");
558 ctmp.Clear("D");
559 }
560 delete htempC;
561 }
562 }
563 psname.str("");
564 psname << Form("cdcdedx_1dcell_dedxhit%s.pdf]", m_suffix.data());
565 ctmp.Print(psname.str().c_str());
566}
567
568//--------------------------------------------------
569void CDCDedx1DCellAlgorithm::plotLayerDist(std::array<TH1D*, 3>& hdedxlay)
570{
571
572 TCanvas cdedxlayer("layerdedxhit", "Inner and Outer Layer dedxhit dist", 900, 400);
573 cdedxlayer.Divide(3, 1);
574
575 for (int il = 0; il < m_kNGroups; il++) {
576 int minlay = 0, maxlay = 0;
577 if (isFixTrunc) {
578 minlay = std::stoi(hdedxlay[il]->GetXaxis()->GetTitle());
579 maxlay = std::stoi(hdedxlay[il]->GetYaxis()->GetTitle());
580 double lowedge = hdedxlay[il]->GetXaxis()->GetBinLowEdge(minlay);
581 double upedge = hdedxlay[il]->GetXaxis()->GetBinUpEdge(maxlay);
582 hdedxlay[il]->SetTitle(Form("%s, trunc #rightarrow: %0.02f - %0.02f;dedxhit;entries", hdedxlay[il]->GetTitle(), lowedge, upedge));
583 }
584
585 cdedxlayer.cd(il + 1);
586 hdedxlay[il]->SetFillColor(kYellow);
587 hdedxlay[il]->Draw("histo");
588
589 if (isFixTrunc) {
590 TH1D* hdedxlayC = (TH1D*)hdedxlay[il]->Clone(Form("hdedxlayC%d", il));
591 hdedxlayC->GetXaxis()->SetRange(minlay, maxlay);
592 hdedxlayC->SetFillColor(kAzure + 1);
593 hdedxlayC->Draw("same histo");
594 }
595 }
596
597 cdedxlayer.SaveAs(Form("cdcdedx_1dcell_dedxlayer_%s.pdf", m_suffix.data()));
598 cdedxlayer.SaveAs(Form("cdcdedx_1dcell_dedxlayer_%s.root", m_suffix.data()));
599}
600
601//--------------------------------------------------
602void CDCDedx1DCellAlgorithm::plotQaPars(std::array<TH1D*, 3>& hentalay, TH2D* hptcosth)
603{
604
605 TCanvas ceadist("ceadist", "Enta distributions", 1600, 800);
606 ceadist.Divide(3, 1);
607
608 for (int il = 0; il < m_kNGroups; il++) {
609
610 ceadist.cd(il + 1);
611 gPad->SetLogy();
612 hentalay[il]->SetFillColor(kYellow);
613 hentalay[il]->Draw("hist");
614 }
615
616 TCanvas cptcos("cptcos", "pt vs costh dist.", 400, 400);
617 cptcos.cd();
618 hptcosth->Draw("colz");
619
620 cptcos.SaveAs(Form("cdcdedx_ptcosth_%s.pdf", m_suffix.data()));
621 ceadist.SaveAs(Form("cdcdedx_1dcell_enta%s.pdf", m_suffix.data()));
622 ceadist.SaveAs(Form("cdcdedx_1dcell_enta%s.root", m_suffix.data()));
623}
624
625//--------------------------------------------------
626void CDCDedx1DCellAlgorithm::plotRelConst(std::vector<double>tempconst, std::vector<double>layerconst, int il)
627{
628
629 TH1D* hconst, *hconstvar;
630
631 Double_t* nvarBins;
632 nvarBins = &m_binValue[il][0];
633
634 hconst = new TH1D(Form("hconst%s", m_label[il].data()), "", m_eaBin, m_eaMin, m_eaMax);
635 std::string title = Form("calibration const dist: %s: (%s); entaRS (#alpha); entries", m_label[il].data(), m_runExp.data());
636 hconst->SetTitle(Form("%s", title.data()));
637
638 hconstvar = new TH1D(Form("hconstvar%s", m_label[il].data()), "", m_eaBinLocal[il], nvarBins);
639 title = Form("calibration const dist (var bins): %s: (%s); entaRS (#alpha);entries", m_label[il].data(), m_runExp.data());
640 hconstvar->SetTitle(Form("%s", title.data()));
641
642 if (isVarBins) {
643 for (int iea = 0; iea < m_eaBinLocal[il]; iea++)
644 hconstvar->SetBinContent(iea + 1, tempconst.at(iea));
645 }
646
647 for (int jea = 0; jea < m_eaBin; jea++) hconst->SetBinContent(jea + 1, layerconst.at(jea));
648
649 gStyle->SetOptStat("ne");
650 TCanvas cconst("cconst", "calibration Constants", 800, 400);
651 if (isVarBins) {
652 cconst.Divide(2, 1);
653 cconst.SetWindowSize(1000, 800);
654 }
655
656 cconst.cd(1);
657 hconst->SetFillColor(kYellow);
658 hconst->Draw("histo");
659 if (isVarBins) {
660 cconst.cd(2);
661 hconstvar->SetFillColor(kBlue);
662 hconstvar->Draw("hist");
663 }
664 cconst.SaveAs(Form("cdcdedx_1dcell_relconst%s_%s.pdf", m_label[il].data(), m_suffix.data()));
665 cconst.SaveAs(Form("cdcdedx_1dcell_relconst%s_%s.root", m_label[il].data(), m_suffix.data()));
666
667 delete hconst;
668 delete hconstvar;
669}
670
671//--------------------------------------------------
673{
674
675 //Draw New/Old final constants
676 TH1D* hnewconst[m_kNGroups], *holdconst[m_kNGroups];
677 double min[m_kNGroups], max[m_kNGroups];
678
679 for (unsigned int il = 0; il < m_kNGroups; il++) {
680 unsigned int nbins = m_DBOneDCell->getNBins(getRepresentativeLayer(il));
681
682 std::string title = Form("final calibration const dist (%s): %s; entaRS (#alpha); entries", m_label[il].data(), m_runExp.data());
683 hnewconst[il] = new TH1D(Form("hnewconst_%s", m_label[il].data()), "", m_eaBin, m_eaMin, m_eaMax);
684 hnewconst[il]->SetTitle(Form("%s", title.data()));
685
686 title = Form("old calibration const dist (%s): %s; entaRS (#alpha); entries", m_label[il].data(), m_runExp.data());
687 holdconst[il] = new TH1D(Form("holdconst_%s", m_label[il].data()), "", m_eaBin, m_eaMin, m_eaMax);
688 holdconst[il]->SetTitle(Form("%s", title.data()));
689
690 for (unsigned int iea = 0; iea < nbins; iea++) {
691 double prev = m_DBOneDCell->getMean(getRepresentativeLayer(il), iea);
692 holdconst[il]->SetBinContent(iea + 1, prev);
693 hnewconst[il]->SetBinContent(iea + 1, m_onedcors[il][iea]);
694 }
695 min[il] = hnewconst[il]->GetMinimum();
696 max[il] = hnewconst[il]->GetMaximum();
697 }
698
699 //Ploting final constants
700 double globalMin = min[0];
701 double globalMax = max[0];
702 for (int il = 1; il < m_kNGroups; il++) {
703 if (min[il] < globalMin) globalMin = min[il];
704 if (max[il] > globalMax) globalMax = max[il];
705 }
706
707 gStyle->SetOptStat("ne");
708 TCanvas cfconst("cfconst", "Final calibration constants", 1600, 600);
709 cfconst.Divide(3, 1);
710
711 for (int il = 0; il < m_kNGroups; il++) {
712 cfconst.cd(il + 1);
713 hnewconst[il]->GetYaxis()->SetRangeUser(globalMin * 0.95, globalMax * 1.05);
714 hnewconst[il]->SetLineColor(kBlack);
715 hnewconst[il]->Draw("histo");
716 holdconst[il]->SetLineColor(kRed);
717 holdconst[il]->Draw("histo same");
718
719 auto legend = new TLegend(0.4, 0.75, 0.56, 0.85);
720 legend->AddEntry(holdconst[il], "Old", "lep");
721 legend->AddEntry(hnewconst[il], "New", "lep");
722 legend->Draw();
723 }
724
725 cfconst.SaveAs(Form("cdcdedx_1dcell_fconsts%s.pdf", m_suffix.data()));
726 cfconst.SaveAs(Form("cdcdedx_1dcell_fconsts%s.root", m_suffix.data()));
727
728 for (int il = 0; il < m_kNGroups; il++) {
729 delete hnewconst[il];
730 delete holdconst[il];
731 }
732}
733
734//------------------------------------
736{
737
738 TCanvas cstats("cstats", "cstats", 1000, 500);
739 cstats.SetBatch(kTRUE);
740 cstats.Divide(2, 1);
741
742 cstats.cd(1);
743 auto hestats = getObjectPtr<TH1I>("hestats");
744 if (hestats) {
745 hestats->SetName(Form("hestats_%s", m_suffix.data()));
746 hestats->SetStats(0);
747 hestats->DrawCopy("");
748 }
749
750 cstats.cd(2);
751 auto htstats = getObjectPtr<TH1I>("htstats");
752 if (htstats) {
753 htstats->SetName(Form("htstats_%s", m_suffix.data()));
754 htstats->SetStats(0);
755 htstats->DrawCopy("");
756 }
757 cstats.Print(Form("cdcdedx_1dcell_stats_%s.pdf", m_suffix.data()));
758}
void plotLayerDist(std::array< TH1D *, m_kNGroups > &hdedxlay)
function to draw dedx dist.
CDCDedx1DCellAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
double m_eaMax
upper edge of entrance angle
void plotMergeFactor(std::map< int, std::vector< double > > bounds, const std::array< int, 2 > nDev, std::map< int, std::vector< int > > steps)
function to plot merging factor
double m_truncMax
upper threshold on truncation
int m_binSplit
multiply nbins by this factor in full range
double m_truncMin
lower threshold on truncation
double m_adjustFac
factor with that one what to adjust baseline
void getTruncatedBins(TH1D *hist, int &binlow, int &binhigh)
function to get bins of truncation from histogram
void CreateBinMapping()
class function to create vectors for bin mapping (Var->symm)
double m_chargeType
charge type for baseline adj
void getExpRunInfo()
function to extract calibration run/exp
unsigned int getRepresentativeLayer(unsigned int il) const
Representative CDC layer for each SL group (used to access group-wise constants): SL0 => 1,...
static constexpr int m_kNGroups
SL grouping: inner (SL0), middle (SL1), outer (SL2–8)
DBObjPtr< CDCDedx1DCell > m_DBOneDCell
One cell correction DB object.
double m_cosMax
a limit on cos theta
bool isPrintLog
print more debug information
std::string m_suffix
add suffix to all plot name
int m_eaB
reset # of bins for entrance angle for each experiment
double getTruncationMean(TH1D *hist, int binlow, int binhigh)
function to get truncated mean
double m_ptMax
a limit on transverse momentum
void plotConstants()
function to draw the old/new final constants
std::array< std::vector< int >, m_kNGroups > m_binIndex
symm/Var bin numbers
bool isFixTrunc
true = fix window for all out/inner layers
bool isVarBins
true: if variable bin size is requested
void plotQaPars(std::array< TH1D *, m_kNGroups > &hentalay, TH2D *hptcosth)
function to draw pt vs costh and entrance angle distribution for Inner/Outer layer
double m_eaBW
binwdith of entrance angle bin
bool isRotSymm
if rotation symmetry requested
std::string m_runExp
add run and exp to title of plot
void defineHisto(std::array< std::vector< TH1D * >, 3 > &hdedxhit, std::array< TH1D *, 3 > &hdedxlay, std::array< TH1D *, 3 > &hentalay)
function to define histograms
void plotEventStats()
function to draw the stats plots
int rotationalBin(int nbin, int ibin)
class function to set rotation symmetry
virtual EResult calibrate() override
1D cell algorithm
std::array< std::vector< double >, m_kNGroups > m_binValue
enta Var bin values
double m_dedxMax
upper edge of dedxhit
void createPayload()
function to generate final constants
bool isMakePlots
produce plots for status
void plotRelConst(std::vector< double >tempconst, std::vector< double >layerconst, int il)
function to draw symm/Var layer constant
std::vector< std::vector< double > > m_onedcors
final vectors of calibration
bool isMerge
print more debug information
double m_dedxMin
lower edge of dedxhit
std::string m_label[m_kNGroups]
add inner/outer superlayer label
void plotdedxHist(std::array< std::vector< TH1D * >, m_kNGroups > &hdedxhit)
function to draw the dE/dx histogram in enta bins
double m_eaMin
lower edge of entrance angle
dE/dx 1D cell correction 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.