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