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_eaB(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 m_eaBinLocal.clear();
214 for (int il = 0; il < 2; il++) {
215 m_binValue[il].clear();
216 m_binIndex[il].clear();
217 }
218 return c_OK;
219}
220
221//--------------------------------------------------
223{
224
225 int cruns = 0;
226 for (auto expRun : getRunList()) {
227 if (cruns == 0) B2INFO("CDCDedxBadWires: start exp " << expRun.first << " and run " << expRun.second << "");
228 cruns++;
229 }
230
231 const auto erStart = getRunList()[0];
232 int estart = erStart.first;
233 int rstart = erStart.second;
234
235 const auto erEnd = getRunList()[cruns - 1];
236 int eend = erEnd.first;
237 int rend = erEnd.second;
238
239 updateDBObjPtrs(1, rstart, estart);
240
241 m_runExp = Form("Range (%d:%d,%d:%d)", estart, rstart, eend, rend);
242 if (m_suffix.length() > 0) m_suffix = Form("%s_e%d_r%dr%d", m_suffix.data(), estart, rstart, rend);
243 else m_suffix = Form("e%d_r%dr%d", estart, rstart, rend);
244}
245
246//--------------------------------------------------
248{
249
250 std::map<int, std::vector<double>> bounds;
251 std::map<int, std::vector<int>> steps;
252
253 const std::array<int, 2> nDev{8, 4};
254 bounds[0] = {0, 108, 123, 133, 158, 183, 193, 208, 316}; //il boundaries
255 steps[0] = {9, 3, 2, 1, 1, 2, 3, 9}; //il steps
256 bounds[1] = {0, 38, 158, 278, 316}; //OL boundaries
257 steps[1] = {2, 1, 1, 2}; //OL steps
258
259 for (int il = 0; il < 2; il++) {
260
261 for (int ibin = 0; ibin <= nDev[il]; ibin++) bounds[il][ibin] = bounds[il][ibin] * m_binSplit;
262
263 int ieaprime = -1, temp = -99, ibin = 0;
264
265 double pastbin = m_eaMin;
266 m_binValue[il].push_back(pastbin);
267
268 for (int iea = 0; iea < m_eaBin; iea++) {
269
270 if (isVarBins) {
271 if (iea % int(bounds[il][ibin + 1]) == 0 && iea > 0) ibin++;
272 int diff = iea - int(bounds[il][ibin]);
273 if (diff % steps[il][ibin] == 0) ieaprime++;
274 } else ieaprime = iea;
275
276 m_binIndex[il].push_back(ieaprime);
277
278 if (ieaprime != temp) {
279 double binwidth = m_eaBW;
280 if (isVarBins) binwidth = m_eaBW * steps[il][ibin];
281 double binvalue = pastbin + binwidth;
282 pastbin = binvalue;
283 if (std::abs(binvalue) < 1e-5)binvalue = 0;
284 m_binValue[il].push_back(binvalue);
285 }
286 temp = ieaprime;
287 }
288 m_eaBinLocal.push_back(int(m_binValue[il].size()) - 1) ;
289 }
290 if (isMakePlots) plotMergeFactor(bounds, nDev, steps);
291}
292
293//--------------------------------------------------
294void CDCDedx1DCellAlgorithm::defineHisto(std::vector<TH1D*> hdedxhit[2], TH1D* hdedxlay[2], TH1D* hentalay[2])
295{
296 for (int il = 0; il < 2; il++) {
297
298 std::string title = Form("dedxhit dist (%s): %s ; dedxhit;entries", m_label[il].data(), m_runExp.data());
299 hdedxlay[il] = new TH1D(Form("hdedxlay%s", m_label[il].data()), "", m_dedxBin, m_dedxMin, m_dedxMax);
300 hdedxlay[il]->SetTitle(Form("%s", title.data()));
301
302 Double_t* nvarBins;
303 nvarBins = &m_binValue[il][0];
304
305 if (isVarBins) title = Form("entaRS dist (variable bins): %s: (%s); entaRS (#alpha);entries", m_label[il].data(), m_runExp.data());
306 else title = Form("entaRS dist (sym. bins): %s: (%s); entaRS (#alpha);entries", m_label[il].data(), m_runExp.data());
307
308 hentalay[il] = new TH1D(Form("hentalay%s", m_label[il].data()), "", m_eaBinLocal[il], nvarBins);
309 hentalay[il]->SetTitle(Form("%s", title.data()));
310
311 for (int iea = 0; iea < m_eaBinLocal[il]; iea++) {
312
313 double min = m_binValue[il].at(iea);
314 double max = m_binValue[il].at(iea + 1);
315 double width = max - min;
316
317 if (isPrintLog) B2INFO("bin: " << iea << " ], min:" << min << " , max: " << max << " , width: " << width);
318
319 title = Form("%s: entaRS = (%0.03f to %0.03f)", m_label[il].data(), min, max);
320 hdedxhit[il].push_back(new TH1D(Form("hdedxhit_%s_bin%d", m_label[il].data(), iea), "", m_dedxBin, m_dedxMin, m_dedxMax));
321 hdedxhit[il][iea]->SetTitle(Form("%s", title.data()));
322 }
323 }
324}
325
326//--------------------------------------------------
327void CDCDedx1DCellAlgorithm::getTruncatedBins(TH1D* hist, int& binlow, int& binhigh)
328{
329
330 //calculating truncation average
331 double sum = hist->Integral();
332 if (sum <= 0 || hist->GetNbinsX() <= 0) {
333 binlow = 1; binhigh = 1;
334 return ;
335 }
336
337 binlow = 1.0; binhigh = 1.0;
338 double sumPer5 = 0.0, sumPer75 = 0.0;
339 for (int ibin = 1; ibin <= hist->GetNbinsX(); ibin++) {
340 double bcdedx = hist->GetBinContent(ibin);
341 if (sumPer5 <= m_truncMin * sum) {
342 sumPer5 += bcdedx;
343 binlow = ibin;
344 }
345 if (sumPer75 <= m_truncMax * sum) {
346 sumPer75 += bcdedx;
347 binhigh = ibin;
348 }
349 }
350 return;
351}
352
353//--------------------------
354double CDCDedx1DCellAlgorithm::getTruncationMean(TH1D* hist, int binlow, int binhigh)
355{
356
357 //calculating truncation average
358 if (hist->Integral() < 100) return 1.0;
359
360 if (binlow <= 0 || binhigh > hist->GetNbinsX())return 1.0;
361
362 double binweights = 0., sumofbc = 0.;
363 for (int ibin = binlow; ibin <= binhigh; ibin++) {
364 double bcdedx = hist->GetBinContent(ibin);
365 if (bcdedx > 0) {
366 binweights += (bcdedx * hist->GetBinCenter(ibin));
367 sumofbc += bcdedx;
368 }
369 }
370 if (sumofbc > 0) return binweights / sumofbc;
371 else return 1.0;
372}
373
374//--------------------------------------------------
376{
377
378 B2INFO("dE/dx one cell calibration: Generating payloads");
379
380 for (unsigned int il = 0; il < 2; il++) {
381 if (isMerge) {
382 unsigned int nbins = m_DBOneDCell->getNBins(il);
383
384 if (int(nbins) != m_eaBin)
385 B2ERROR("merging failed because of unmatch bins (old " << m_eaBin << " new " << nbins << ")");
386
387 for (unsigned int iea = 0; iea < nbins; iea++) {
388 double prev = m_DBOneDCell->getMean(8 * il + 1, iea);
389 m_onedcors[il][iea] *= prev;
390 // m_onedcors[il][iea] /= 0.98;
391 }
392 }
393
394 if (m_chargeType > 0)
395 for (int ie = 0; ie < m_eaBin / 2; ie++) m_onedcors[il][ie] *= m_adjustFac;
396 if (m_chargeType < 0)
397 for (int ie = m_eaBin / 2; ie < m_eaBin; ie++) m_onedcors[il][ie] *= m_adjustFac;
398
399 }
400 //Saving constants
401 B2INFO("dE/dx Calibration done for CDCDedx1DCell");
402 CDCDedx1DCell* gain = new CDCDedx1DCell(0, m_onedcors);
403 saveCalibration(gain, "CDCDedx1DCell");
404}
405
406//--------------------------------------------------
407void CDCDedx1DCellAlgorithm::plotMergeFactor(std::map<int, std::vector<double>> bounds, const std::array<int, 2> nDev,
408 std::map<int, std::vector<int>> steps)
409{
410
411 TCanvas cmfactor("cmfactor", "Merging factors", 800, 400);
412 cmfactor.Divide(2, 1);
413
414 for (int il = 0; il < 2; il++) {
415 Double_t* nvarBins;
416 nvarBins = &bounds[il][0];
417
418 TH1I* hist = new TH1I(Form("hist_%s", m_label[il].data()), "", nDev[il], nvarBins);
419 hist->SetTitle(Form("Merging factor for %s bins;binindex;merge-factors", m_label[il].data()));
420
421 for (int ibin = 0; ibin < nDev[il]; ibin++) hist->SetBinContent(ibin + 1, steps[il][ibin]);
422
423 cmfactor.cd(il + 1);
424 hist->SetFillColor(kYellow);
425 hist->Draw("hist");
426 delete hist;
427 }
428
429 cmfactor.SaveAs(Form("cdcdedx_1dcell_mergefactor%s.pdf", m_suffix.data()));
430 cmfactor.SaveAs(Form("cdcdedx_1dcell_mergefactor%s.root", m_suffix.data()));
431}
432
433//--------------------------------------------------
434void CDCDedx1DCellAlgorithm::plotdedxHist(std::vector<TH1D*> hdedxhit[2])
435{
436
437 TCanvas ctmp("tmp", "tmp", 1200, 1200);
438 ctmp.Divide(4, 4);
439 std::stringstream psname;
440
441 psname << Form("cdcdedx_1dcell_dedxhit%s.pdf[", m_suffix.data());
442 ctmp.Print(psname.str().c_str());
443 psname.str("");
444 psname << Form("cdcdedx_1dcell_dedxhit%s.pdf", m_suffix.data());
445
446 for (int il = 0; il < 2; il++) {
447
448 for (int jea = 0; jea < m_eaBinLocal[il]; jea++) {
449
450 int minbin = std::stoi(hdedxhit[il][jea]->GetXaxis()->GetTitle());
451 int maxbin = std::stoi(hdedxhit[il][jea]->GetYaxis()->GetTitle());
452
453 ctmp.cd(jea % 16 + 1);
454 hdedxhit[il][jea]->SetFillColor(4 + il);
455
456 hdedxhit[il][jea]->SetTitle(Form("%s;dedxhit;entries", hdedxhit[il][jea]->GetTitle()));
457 hdedxhit[il][jea]->DrawClone("hist");
458 TH1D* htempC = (TH1D*)hdedxhit[il][jea]->Clone(Form("%sc2", hdedxhit[il][jea]->GetName()));
459 htempC->GetXaxis()->SetRange(minbin, maxbin);
460 htempC->SetFillColor(kGray);
461 htempC->DrawClone("same hist");
462
463 if (jea % 16 == 15 || (jea == m_eaBinLocal[il] - 1)) {
464 ctmp.Print(psname.str().c_str());
465 gPad->Clear("D");
466 ctmp.Clear("D");
467 }
468 delete htempC;
469 }
470 }
471 psname.str("");
472 psname << Form("cdcdedx_1dcell_dedxhit%s.pdf]", m_suffix.data());
473 ctmp.Print(psname.str().c_str());
474}
475
476//--------------------------------------------------
478{
479
480 TCanvas cdedxlayer("layerdedxhit", "Inner and Outer Layer dedxhit dist", 900, 400);
481 cdedxlayer.Divide(2, 1);
482
483 for (int il = 0; il < 2; il++) {
484 int minlay = 0, maxlay = 0;
485 if (isFixTrunc) {
486 minlay = std::stoi(hdedxlay[il]->GetXaxis()->GetTitle());
487 maxlay = std::stoi(hdedxlay[il]->GetYaxis()->GetTitle());
488 double lowedge = hdedxlay[il]->GetXaxis()->GetBinLowEdge(minlay);
489 double upedge = hdedxlay[il]->GetXaxis()->GetBinUpEdge(maxlay);
490 hdedxlay[il]->SetTitle(Form("%s, trunc #rightarrow: %0.02f - %0.02f;dedxhit;entries", hdedxlay[il]->GetTitle(), lowedge, upedge));
491 }
492
493 cdedxlayer.cd(il + 1);
494 hdedxlay[il]->SetFillColor(kYellow);
495 hdedxlay[il]->Draw("histo");
496
497 if (isFixTrunc) {
498 TH1D* hdedxlayC = (TH1D*)hdedxlay[il]->Clone(Form("hdedxlayC%d", il));
499 hdedxlayC->GetXaxis()->SetRange(minlay, maxlay);
500 hdedxlayC->SetFillColor(kAzure + 1);
501 hdedxlayC->Draw("same histo");
502 }
503 }
504
505 cdedxlayer.SaveAs(Form("cdcdedx_1dcell_dedxlay%s.pdf", m_suffix.data()));
506 cdedxlayer.SaveAs(Form("cdcdedx_1dcell_dedxlay%s.root", m_suffix.data()));
507}
508
509//--------------------------------------------------
510void CDCDedx1DCellAlgorithm::plotQaPars(TH1D* hentalay[2], TH2D* hptcosth)
511{
512
513 TCanvas ceadist("ceadist", "Enta distributions", 800, 400);
514 ceadist.Divide(2, 1);
515
516 for (int il = 0; il < 2; il++) {
517
518 ceadist.cd(il + 1);
519 gPad->SetLogy();
520 hentalay[il]->SetFillColor(kYellow);
521 hentalay[il]->Draw("hist");
522 }
523
524 TCanvas cptcos("cptcos", "pt vs costh dist.", 400, 400);
525 cptcos.cd();
526 hptcosth->Draw("colz");
527
528 cptcos.SaveAs(Form("cdcdedx_ptcosth_%s.pdf", m_suffix.data()));
529 ceadist.SaveAs(Form("cdcdedx_1dcell_enta%s.pdf", m_suffix.data()));
530 ceadist.SaveAs(Form("cdcdedx_1dcell_enta%s.root", m_suffix.data()));
531}
532
533//--------------------------------------------------
534void CDCDedx1DCellAlgorithm::plotRelConst(std::vector<double>tempconst, std::vector<double>layerconst, int il)
535{
536
537 TH1D* hconst, *hconstvar;
538
539 Double_t* nvarBins;
540 nvarBins = &m_binValue[il][0];
541
542 hconst = new TH1D(Form("hconst%s", m_label[il].data()), "", m_eaBin, m_eaMin, m_eaMax);
543 std::string title = Form("calibration const dist: %s: (%s); entaRS (#alpha); entries", m_label[il].data(), m_runExp.data());
544 hconst->SetTitle(Form("%s", title.data()));
545
546 hconstvar = new TH1D(Form("hconstvar%s", m_label[il].data()), "", m_eaBinLocal[il], nvarBins);
547 title = Form("calibration const dist (var bins): %s: (%s); entaRS (#alpha);entries", m_label[il].data(), m_runExp.data());
548 hconstvar->SetTitle(Form("%s", title.data()));
549
550 if (isVarBins) {
551 for (int iea = 0; iea < m_eaBinLocal[il]; iea++)
552 hconstvar->SetBinContent(iea + 1, tempconst.at(iea));
553 }
554
555 for (int jea = 0; jea < m_eaBin; jea++) hconst->SetBinContent(jea + 1, layerconst.at(jea));
556
557 gStyle->SetOptStat("ne");
558 TCanvas cconst("cconst", "Calirbation Constants", 800, 400);
559 if (isVarBins) {
560 cconst.Divide(2, 1);
561 cconst.SetWindowSize(1000, 800);
562 }
563
564 cconst.cd(1);
565 hconst->SetFillColor(kYellow);
566 hconst->Draw("histo");
567 if (isVarBins) {
568 cconst.cd(2);
569 hconstvar->SetFillColor(kBlue);
570 hconstvar->Draw("hist");
571 }
572 cconst.SaveAs(Form("cdcdedx_1dcell_relconst%s_%s.pdf", m_label[il].data(), m_suffix.data()));
573 cconst.SaveAs(Form("cdcdedx_1dcell_relconst%s_%s.root", m_label[il].data(), m_suffix.data()));
574
575 delete hconst;
576 delete hconstvar;
577}
578
579//--------------------------------------------------
581{
582
583 //Draw New/Old final constants
584 TH1D* hnewconst[2], *holdconst[2];
585 double min[2], max[2];
586
587 for (unsigned int il = 0; il < 2; il++) {
588 unsigned int nbins = m_DBOneDCell->getNBins(il);
589
590 std::string title = Form("final calibration const dist (%s): %s; entaRS (#alpha); entries", m_label[il].data(), m_runExp.data());
591 hnewconst[il] = new TH1D(Form("hnewconst_%s", m_label[il].data()), "", m_eaBin, m_eaMin, m_eaMax);
592 hnewconst[il]->SetTitle(Form("%s", title.data()));
593
594 title = Form("old calibration const dist (%s): %s; entaRS (#alpha); entries", m_label[il].data(), m_runExp.data());
595 holdconst[il] = new TH1D(Form("holdconst_%s", m_label[il].data()), "", m_eaBin, m_eaMin, m_eaMax);
596 holdconst[il]->SetTitle(Form("%s", title.data()));
597
598 for (unsigned int iea = 0; iea < nbins; iea++) {
599 double prev = m_DBOneDCell->getMean(8 * il + 1, iea);
600 holdconst[il]->SetBinContent(iea + 1, prev);
601 hnewconst[il]->SetBinContent(iea + 1, m_onedcors[il][iea]);
602 }
603 min[il] = hnewconst[il]->GetMinimum();
604 max[il] = hnewconst[il]->GetMaximum();
605 }
606
607 //Plotting final constants
608 if (max[1] < max[0])max[1] = max[0];
609 if (min[1] > min[0])min[1] = min[0];
610
611 gStyle->SetOptStat("ne");
612 TCanvas cfconst("cfconst", "Final calirbation constants", 800, 400);
613 cfconst.Divide(2, 1);
614
615 for (int il = 0; il < 2; il++) {
616 cfconst.cd(il + 1);
617 hnewconst[il]->GetYaxis()->SetRangeUser(min[1] * 0.95, max[1] * 1.05);
618 hnewconst[il]->SetLineColor(kBlack);
619 hnewconst[il]->Draw("histo");
620 holdconst[il]->SetLineColor(kRed);
621 holdconst[il]->Draw("histo same");
622
623 auto legend = new TLegend(0.4, 0.75, 0.56, 0.85);
624 legend->AddEntry(holdconst[il], "Old", "lep");
625 legend->AddEntry(hnewconst[il], "New", "lep");
626 legend->Draw();
627 }
628
629 cfconst.SaveAs(Form("cdcdedx_1dcell_fconsts%s.pdf", m_suffix.data()));
630 cfconst.SaveAs(Form("cdcdedx_1dcell_fconsts%s.root", m_suffix.data()));
631
632 for (int il = 0; il < 2; il++) {
633 delete hnewconst[il];
634 delete holdconst[il];
635 }
636}
637
638//------------------------------------
640{
641
642 TCanvas cstats("cstats", "cstats", 1000, 500);
643 cstats.SetBatch(kTRUE);
644 cstats.Divide(2, 1);
645
646 cstats.cd(1);
647 auto hestats = getObjectPtr<TH1I>("hestats");
648 if (hestats) {
649 hestats->SetName(Form("hestats_%s", m_suffix.data()));
650 hestats->SetStats(0);
651 hestats->DrawCopy("");
652 }
653
654 cstats.cd(2);
655 auto htstats = getObjectPtr<TH1I>("htstats");
656 if (htstats) {
657 hestats->SetName(Form("htstats_%s", m_suffix.data()));
658 htstats->DrawCopy("");
659 hestats->SetStats(0);
660 }
661 cstats.Print(Form("cdcdedx_1dcell_stats_%s.pdf", m_suffix.data()));
662}
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 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
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
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
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 entrance 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 entrance 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.