Belle II Software development
CDCDedxCosineAlgorithm.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/CDCDedxCosineAlgorithm.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//-----------------------------------------------------------------
20// Implementation
21//-----------------------------------------------------------------
23 CalibrationAlgorithm("CDCDedxElectronCollector"),
24 isMethodSep(true),
25 isMakePlots(true),
26 isMergePayload(true),
27 m_sigLim(2.5),
28 m_cosBin(100),
29 m_cosMin(-1.0),
30 m_cosMax(1.0),
31 m_dedxBin(250),
32 m_dedxMin(0.0),
33 m_dedxMax(5.0),
34 m_suffix("")
35
36{
37 // Set module properties
38 setDescription("A calibration algorithm for CDC dE/dx electron cos(theta) dependence");
39
40}
41
42//-----------------------------------------------------------------
43// Run the calibration
44//-----------------------------------------------------------------
46{
47
49
50 if (!m_DBCosineCor.isValid())
51 B2FATAL("There is no valid previous payload for CDCDedxCosineCor");
52
53 B2INFO("Preparing dE/dx calibration for CDC dE/dx electron saturation");
54
55 // Get data objects
56 auto ttree = getObjectPtr<TTree>("tree");
57 if (!ttree) {
58 B2ERROR("Input tree 'tree' not found");
59 return c_Failure;
60 }
61 if (ttree->GetEntries() < 100)return c_NotEnoughData;
62
63 double dedx, costh; int charge;
64 ttree->SetBranchAddress("dedx", &dedx);
65 ttree->SetBranchAddress("costh", &costh);
66 ttree->SetBranchAddress("charge", &charge);
67
68
69 // make histograms to store dE/dx values in bins of cos(theta)
70 // bin size can be arbitrary, but for now just make uniform bins
71 std::vector<TH1D*> hDedxCos_neg, hDedxCos_pos, hDedxCos_all;
72
73 const double binW = (m_cosMax - m_cosMin) / m_cosBin;
74
75 defineHisto(hDedxCos_neg, "neg", "e-");
76 defineHisto(hDedxCos_pos, "pos", "e+");
77 defineHisto(hDedxCos_all, "all", "e-,e+");
78
79 // fill histograms, bin size may be arbitrary
80 TH1D* hCosth_neg = defineCosthHist("neg");
81 TH1D* hCosth_pos = defineCosthHist("pos");
82 TH1D* hCosth_all = defineCosthHist("all");
83
84 for (int i = 0; i < ttree->GetEntries(); ++i) {
85
86 ttree->GetEvent(i);
87
88 //if track is a junk
89 if (dedx <= 0 || charge == 0) continue;
90
91 //if track is in CDC accpetance (though it is inbuilt in collector module)
92 if (costh < TMath::Cos(150 * TMath::DegToRad()) || costh > TMath::Cos(17 * TMath::DegToRad())) continue;
93
94 int bin = int((costh - m_cosMin) / binW);
95 if (bin < 0 || bin >= static_cast<int>(m_cosBin)) continue;
96
97 if (isMethodSep) {
98 if (charge < 0) {
99 hCosth_neg->Fill(costh);
100 hDedxCos_neg[bin]->Fill(dedx);
101 } else if (charge > 0) {
102 hCosth_pos->Fill(costh);
103 hDedxCos_pos[bin]->Fill(dedx);
104 }
105 } else {
106 hCosth_all->Fill(costh);
107 hDedxCos_all[bin]->Fill(dedx);
108 }
109 }
110
111 // fit histograms to get gains in bins of cos(theta)
112 std::vector<double> cosine;
113
114 std::vector<std::vector<double>> dedxAll(4);
115 std::vector<std::vector<double>> dedxNeg(4);
116 std::vector<std::vector<double>> dedxPos(4);
117
118 for (unsigned int i = 0; i < m_cosBin; ++i) {
119
120 double meanDedx = 1.0; //This is what we need for calibration
121 double meanDedxErr = 0.0;
122
123 if (!isMethodSep) {
124 FitValues fitAll = fitHistogram(hDedxCos_all[i]);
125
126 meanDedx = fitAll.mean;
127
128 dedxAll[0].push_back(fitAll.mean);
129 dedxAll[1].push_back(fitAll.meanErr);
130 dedxAll[2].push_back(fitAll.sigma);
131 dedxAll[3].push_back(fitAll.sigmaErr);
132
133 } else {
134
135 //Fit electron dE/dx in cos bins
136 FitValues fitNeg = fitHistogram(hDedxCos_neg[i]);
137
138 //Fit positron dE/dx in cos bins
139 FitValues fitPos = fitHistogram(hDedxCos_pos[i]);
140
141 if (fitPos.status != "FitOK" && fitNeg.status == "FitOK") {
142 fitPos.mean = fitNeg.mean;
143 hDedxCos_pos[i]->SetTitle(Form("%s, mean (manual) = elec left", hDedxCos_pos[i]->GetTitle()));
144 } else if (fitNeg.status != "FitOK" && fitPos.status == "FitOK") {
145 fitNeg.mean = fitPos.mean;
146 hDedxCos_neg[i]->SetTitle(Form("%s, mean (manual) = posi right", hDedxCos_neg[i]->GetTitle()));
147 } else if (fitNeg.status != "FitOK" && fitPos.status != "FitOK") {
148 fitNeg.mean = 1.0;
149 fitPos.mean = 1.0;
150 }
151
152 dedxNeg[0].push_back(fitNeg.mean);
153 dedxNeg[1].push_back(fitNeg.meanErr);
154 dedxNeg[2].push_back(fitNeg.sigma);
155 dedxNeg[3].push_back(fitNeg.sigmaErr);
156
157 dedxPos[0].push_back(fitPos.mean);
158 dedxPos[1].push_back(fitPos.meanErr);
159 dedxPos[2].push_back(fitPos.sigma);
160 dedxPos[3].push_back(fitPos.sigmaErr);
161
162 meanDedx = 0.5 * (fitNeg.mean + fitPos.mean);
163 if (meanDedx <= 0.0) meanDedx = 1.0;
164
165 meanDedxErr = 0.5 * TMath::Sqrt(fitNeg.meanErr * fitNeg.meanErr +
166 fitPos.meanErr * fitPos.meanErr);
167
168 dedxAll[0].push_back(meanDedx);
169 dedxAll[1].push_back(meanDedxErr);
170
171 }
172
173 cosine.push_back(meanDedx);
174 }
175
176 createPayload(cosine);
177
178
179 if (isMakePlots) {
180
181 //1. dE/dx dist. for cosine bins
182 plotdedxHist(hDedxCos_all, hDedxCos_neg, hDedxCos_pos);
183
184 //4. costh distribution
185 plotCosThetaDist(hCosth_all, hCosth_pos, hCosth_neg);
186
187 plotFitResults(dedxAll, dedxNeg, dedxPos);
188
189 //7. plot statistics related plots here
191
192 //6. draw the final constants
194 }
195
196 m_suffix.clear();
197
198 return c_OK;
199}
200
201//--------------------------------------------------
203{
204
205 int cruns = 0;
206 for (auto expRun : getRunList()) {
207 if (cruns == 0) B2INFO("CDCDedxCosineCor: start exp " << expRun.first << " and run " << expRun.second << "");
208 cruns++;
209 }
210
211 const auto erStart = getRunList()[0];
212 int estart = erStart.first;
213 int rstart = erStart.second;
214
215 const auto erEnd = getRunList()[cruns - 1];
216 int rend = erEnd.second;
217
218 updateDBObjPtrs(1, rstart, estart);
219
220 if (m_suffix.length() > 0) m_suffix = Form("%s_e%d_r%dr%d", m_suffix.data(), estart, rstart, rend);
221 else m_suffix = Form("e%d_r%dr%d", estart, rstart, rend);
222}
223
224//--------------------------------------------------
225void CDCDedxCosineAlgorithm::defineHisto(std::vector<TH1D*>& hdedx, const std::string& tag,
226 const std::string& chargeLabel)
227{
228
229 const double binW = (m_cosMax - m_cosMin) / m_cosBin;
230
231 hdedx.reserve(m_cosBin);
232
233 for (unsigned int i = 0; i < m_cosBin; ++i) {
234 double coslow = i * binW + m_cosMin;
235 double coshigh = coslow + binW;
236
237 hdedx.push_back(new TH1D(Form("hDedxCos_%s_bin%d_%s", tag.c_str(), i, m_suffix.data()), "", m_dedxBin, m_dedxMin, m_dedxMax));
238
239 hdedx[i]->SetTitle(Form("dE/dx dist (%s) in costh (%0.02f, %0.02f);dE/dx (no had sat, for %s);Entries", chargeLabel.c_str(), coslow,
240 coshigh, chargeLabel.c_str()));
241
242 }
243}
244
245//--------------------------------------------------
246TH1D* CDCDedxCosineAlgorithm::defineCosthHist(const std::string& tag)
247{
248
249 TH1D* hist = new TH1D(Form("hCosth_%s_%s", tag.c_str(), m_suffix.data()), " ", m_cosBin, m_cosMin, m_cosMax);
250 hist->SetTitle("cos(#theta) dist (e- and e+); cos(#theta); Entries");
251
252 return hist;
253}
254
255//----------------------------------------
257{
258 FitValues fitValues;
259
260 fitGaussianWithRange(hist, fitValues.status);
261
262 hist->SetFillColorAlpha(kAzure + 1, 0.30);
263
264 if (fitValues.status == "FitOK") {
265 TF1* fitFunc = hist->GetFunction("gaus");
266 if (fitFunc) {
267 fitValues.mean = fitFunc->GetParameter(1);
268 fitValues.meanErr = fitFunc->GetParError(1);
269 fitValues.sigma = fitFunc->GetParameter(2);
270 fitValues.sigmaErr = fitFunc->GetParError(2);
271
272 std::string fitSummary = Form("#mu_{fit}: %0.03f #pm %0.03f, #sigma_{fit}: %0.03f",
273 fitValues.mean, fitValues.meanErr, fitValues.sigma);
274
275 hist->SetTitle(Form("%s, %s", hist->GetTitle(), fitSummary.data()));
276 }
277 }
278
279 return fitValues;
280}
281
282//--------------------------------------------------
283void CDCDedxCosineAlgorithm::fitGaussianWithRange(TH1D*& temphist, TString& status)
284{
285 if (temphist->Integral() < 2000) { //atleast 1k bhabha events
286 B2INFO(Form("\tThis hist (%s) have insufficient entries to perform fit (%0.03f)", temphist->GetName(), temphist->Integral()));
287 status = "LowStats";
288 return;
289 } else {
290 temphist->GetXaxis()->SetRange(temphist->FindFirstBinAbove(0, 1), temphist->FindLastBinAbove(0, 1));
291 int fs = temphist->Fit("gaus", "QR");
292 if (fs != 0) {
293 B2INFO(Form("\tFit (round 1) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
294 status = "FitFailed";
295 return;
296 } else {
297 double meanDedx = temphist->GetFunction("gaus")->GetParameter(1);
298 double width = temphist->GetFunction("gaus")->GetParameter(2);
299 temphist->GetXaxis()->SetRangeUser(meanDedx - 5.0 * width, meanDedx + 5.0 * width);
300 fs = temphist->Fit("gaus", "QR", "", meanDedx - m_sigLim * width, meanDedx + m_sigLim * width);
301 if (fs != 0) {
302 B2INFO(Form("\tFit (round 2) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
303 status = "FitFailed";
304 return;
305 } else {
306 temphist->GetXaxis()->SetRangeUser(meanDedx - 5.0 * width, meanDedx + 5.0 * width);
307 B2INFO(Form("\tFit for hist (%s) sucessfull (status = %d)", temphist->GetName(), fs));
308 status = "FitOK";
309 }
310 }
311 }
312}
313
314//--------------------------------------------------
315void CDCDedxCosineAlgorithm::createPayload(std::vector<double> cosine)
316{
317 m_coscors.resize(m_kNGroups);
318
319 for (unsigned int il = 0; il < m_kNGroups; il++) {
320
321 unsigned int nbins = m_DBCosineCor->getSize(getRepresentativeLayer(il));
322
323 if (nbins != m_cosBin)
324 B2ERROR("merging failed because of unmatch bins (old "
325 << nbins << " new " << m_cosBin << ")");
326
327 m_coscors[il].reserve(nbins);
328
329 for (unsigned int ibin = 0; ibin < nbins; ibin++) {
330
331 double value = cosine[ibin];
332
333 if (isMergePayload) {
334 double prev = m_DBCosineCor->getMean(getRepresentativeLayer(il), ibin);
335
336 value *= prev;
337
338 B2INFO("Cosine Corr for " << m_label[il]
339 << " Bin # " << ibin
340 << ", Previous = " << prev
341 << ", Relative = " << cosine[ibin]
342 << ", Merged = " << value);
343 }
344
345 m_coscors[il].push_back(value);
346 }
347 }
348
349 //Saving constants
350 B2INFO("dE/dx calibration done for CDC dE/dx electron saturation");
351
352 std::vector<unsigned int> layerToGroup(56);
353
354 for (unsigned int layer = 0; layer < 56; layer++) {
355 if (layer < 8) layerToGroup[layer] = 0; // SL0
356 else if (layer < 14) layerToGroup[layer] = 1; // SL1
357 else layerToGroup[layer] = 2; // SL2-8
358 }
359
360 CDCDedxCosineCor* gain = new CDCDedxCosineCor(m_coscors, layerToGroup);
361 saveCalibration(gain, "CDCDedxCosineCor");
362}
363
364//--------------------------------------------------
365void CDCDedxCosineAlgorithm::plotdedxHist(std::vector<TH1D*>& hDedxCos_all,
366 std::vector<TH1D*>& hDedxCos_neg,
367 std::vector<TH1D*>& hDedxCos_pos)
368{
369
370 TCanvas ctmp("tmp", "tmp", 1200, 1200);
371 int nx = isMethodSep ? 2 : 2;
372 int ny = isMethodSep ? 1 : 2;
373 unsigned int nPads = nx * ny;
374 if (isMethodSep) ctmp.SetCanvasSize(1200, 600);
375 ctmp.Divide(nx, ny);
376 std::stringstream psname;
377
378 psname << Form("cdcdedx_coscorr_dedx_%s.pdf[", m_suffix.data());
379 ctmp.Print(psname.str().c_str());
380 psname.str("");
381 psname << Form("cdcdedx_coscorr_dedx_%s.pdf", m_suffix.data());
382
383 for (unsigned int ic = 0; ic < m_cosBin; ic++) {
384 if (!isMethodSep) {
385 ctmp.cd(ic % nPads + 1);
386 hDedxCos_all[ic]->SetStats(0);
387 hDedxCos_all[ic]->SetFillColorAlpha(kYellow, 0.25);
388 hDedxCos_all[ic]->DrawCopy();
389
390 if (ic % nPads == nPads - 1 || ic == m_cosBin - 1) {
391 ctmp.Print(psname.str().c_str());
392 ctmp.Clear();
393 ctmp.Divide(nx, ny);
394 }
395 } else {
396
397 // left: electron
398 ctmp.cd(1);
399 hDedxCos_neg[ic]->SetFillColorAlpha(kRed, 0.25);
400 hDedxCos_neg[ic]->DrawCopy();
401
402
403 // right: positron
404 ctmp.cd(2);
405 hDedxCos_pos[ic]->SetFillColorAlpha(kBlue, 0.25);
406 hDedxCos_pos[ic]->DrawCopy();
407
408 ctmp.Print(psname.str().c_str());
409 ctmp.Clear();
410 ctmp.Divide(nx, ny);
411
412 }
413
414 }
415 psname.str("");
416 psname << Form("cdcdedx_coscorr_dedx_%s.pdf]", m_suffix.data());
417 ctmp.Print(psname.str().c_str());
418}
419
420//--------------------------------------------------
421void CDCDedxCosineAlgorithm::plotCosThetaDist(TH1D* hCosth_all, TH1D* hCosth_pos, TH1D* hCosth_neg)
422{
423
424 TCanvas ceadist("ceadist", "Cosine distributions", 800, 600);
425 ceadist.cd();
426
427
428 // If method separation, overlay pos/neg
429 if (isMethodSep) {
430 TLegend* leg = new TLegend(0.6, 0.7, 0.8, 0.9);
431
432 if (hCosth_neg) {
433 hCosth_neg->SetLineColor(kRed);
434 hCosth_neg->SetFillColorAlpha(kYellow, 0.55);
435 hCosth_neg->SetStats(0);
436 hCosth_neg->Draw("hist");
437 leg->AddEntry(hCosth_neg, "neg", "f");
438 }
439 if (hCosth_pos) {
440 hCosth_pos->SetLineColor(kBlue);
441 hCosth_pos->SetFillColorAlpha(kGray, 0.35);
442 hCosth_pos->SetStats(0);
443 hCosth_pos->Draw("hist same");
444 leg->AddEntry(hCosth_pos, "pos", "f");
445 }
446 leg->Draw();
447
448 } else {
449 // Always draw ALL first
450 if (hCosth_all) {
451 hCosth_all->SetFillColorAlpha(kGray, 0.25);
452 hCosth_all->SetLineColor(kGray);
453 hCosth_all->SetStats(0);
454 hCosth_all->Draw("hist");
455 }
456 }
457
458 ceadist.SaveAs(Form("cdcdedx_coscorr_cosine_%s.pdf", m_suffix.data()));
459 ceadist.SaveAs(Form("cdcdedx_coscorr_cosine_%s.root", m_suffix.data()));
460}
461
462
463//--------------------------------------------------
464void CDCDedxCosineAlgorithm::plotFitResults(const std::vector<std::vector<double>>& dedxAll,
465 const std::vector<std::vector<double>>& dedxNeg,
466 const std::vector<std::vector<double>>& dedxPos)
467{
468 // Fill histograms
469
470 TH1D* hMean_all = new TH1D("hMean_all", "mean vs cos#theta;cos#theta;mean", m_cosBin, m_cosMin, m_cosMax);
471
472 TH1D* hMean_el = new TH1D("hMean_el", "mean (e-);cos#theta;mean", m_cosBin, m_cosMin, m_cosMax);
473
474 TH1D* hMean_po = new TH1D("hMean_po", "mean (e+);cos#theta;mean", m_cosBin, m_cosMin, m_cosMax);
475
476
477 TH1D* hSig_all = new TH1D("hSig_all", "sigma vs cos#theta;cos#theta;#sigma", m_cosBin, m_cosMin, m_cosMax);
478
479 TH1D* hSig_el = new TH1D("hSig_el", "sigma (e-);cos#theta;#sigma", m_cosBin, m_cosMin, m_cosMax);
480
481 TH1D* hSig_po = new TH1D("hSig_po", "sigma (e+);cos#theta;#sigma", m_cosBin, m_cosMin, m_cosMax);
482
483
484 for (unsigned int i = 0; i < m_cosBin; i++) {
485 hMean_all->SetBinContent(i + 1, dedxAll[0][i]);
486 hMean_all->SetBinError(i + 1, dedxAll[1][i]);
487
488 if (isMethodSep) {
489 hMean_el->SetBinContent(i + 1, dedxNeg[0][i]);
490 hMean_el->SetBinError(i + 1, dedxNeg[1][i]);
491 hSig_el->SetBinContent(i + 1, dedxNeg[2][i]);
492 hSig_el->SetBinError(i + 1, dedxNeg[3][i]);
493
494 hMean_po->SetBinContent(i + 1, dedxPos[0][i]);
495 hMean_po->SetBinError(i + 1, dedxPos[1][i]);
496 hSig_po->SetBinContent(i + 1, dedxPos[2][i]);
497 hSig_po->SetBinError(i + 1, dedxPos[3][i]);
498 } else {
499 hSig_all->SetBinContent(i + 1, dedxAll[2][i]);
500 hSig_all->SetBinError(i + 1, dedxAll[3][i]);
501 }
502
503 }
504
505 TCanvas* ctmp = new TCanvas("c_fit", "Mean & Sigma", 1000, 500);
506 ctmp->Divide(2, 1);
507 ctmp->cd(1);
508 gPad->SetGridy(1);
509
510 setHist(hMean_all, kBlack, "dedx rel(#mu_{fit}) for e- and e+ combined", 0.97, 1.04);
511
512 if (isMethodSep) {
513
514 setHist(hMean_el, kRed, "comparison of dedx #mu_{fit}^{rel}", 0.96, 1.04);
515 setHist(hMean_po, kBlue, "", 0.96, 1.04);
516
517 hMean_el->Draw("E1");
518 hMean_po->Draw("E1 same");
519 hMean_all->Draw("E1 same");
520
521 } else {
522 hMean_all->Draw("E1");
523 }
524
525 ctmp->cd(2);
526 gPad->SetGridy(1);
527
528 setHist(hSig_all, kBlack, "dedx rel(#sigma_{fit}) for e- and e+ combined", 0.05, 0.23);
529
530 if (isMethodSep) {
531
532 setHist(hSig_el, kRed, "comparison of dedx #sigma_{fit}^{rel}", 0.05, 0.23, 24);
533 setHist(hSig_po, kBlue, "", 0.05, 0.23, 25);
534
535 hSig_el->Draw("E1");
536 hSig_po->Draw("E1 same");
537
538 } else {
539 hSig_all->Draw("E1");
540 }
541 B2INFO("Plotting finished ");
542
543
544 ctmp->SaveAs(Form("cdcdedx_coscorr_fit_%s.pdf", m_suffix.data()));
545 delete hMean_all;
546 delete hMean_el;
547 delete hMean_po;
548 delete hSig_all;
549 delete hSig_el;
550 delete hSig_po;
551 delete ctmp;
552}
553
554//--------------------------------------------------
556{
557
558 for (int il = 0; il < m_kNGroups; il++) {
559
560 unsigned int nbins = m_DBCosineCor->getSize(getRepresentativeLayer(il));
561
562 // --- Create histograms ---
563 TH1D* hnew = new TH1D(Form("hnew_%s", m_label[il].data()), Form("Final const: %s;cos(#theta);dedx #mu_{fit}", m_label[il].data()),
565 m_cosMax);
566
567 TH1D* hold = new TH1D(Form("hold_%s", m_label[il].data()), Form("Final const: %s;cos(#theta);dedx #mu_{fit}", m_label[il].data()),
569 m_cosMax);
570
571 for (unsigned int iea = 0; iea < nbins; iea++) {
572 double oldv = m_DBCosineCor->getMean(getRepresentativeLayer(il), iea);
573 double newv = m_coscors[il][iea];
574
575 hold->SetBinContent(iea + 1, oldv);
576 hnew->SetBinContent(iea + 1, newv);
577 }
578
579 // --- Ratio ---
580 TH1D* hratio = (TH1D*)hnew->Clone(Form("hratio_%s", m_label[il].data()));
581 hratio->Divide(hold);
582
583 TCanvas c(Form("c_%s", m_label[il].data()), Form("Final constants %s", m_label[il].data()), 1000, 500);
584 c.Divide(2, 1);
585 c.cd(1);
586 gPad->SetGridy(1);
587 gPad->SetGridx(1);
588
589 hnew->SetLineColor(kBlack);
590 hnew->SetStats(0);
591 hold->SetLineColor(kRed);
592 hold->SetStats(0);
593
594 double min = std::min(hnew->GetMinimum(), hold->GetMinimum());
595 double max = std::max(hnew->GetMaximum(), hold->GetMaximum());
596 hnew->GetYaxis()->SetRangeUser(min * 0.95, max * 1.05);
597
598 hnew->Draw("hist");
599 hold->Draw("hist same");
600
601 auto leg = new TLegend(0.6, 0.75, 0.85, 0.88);
602 leg->SetBorderSize(0);
603 leg->SetFillStyle(0);
604 leg->AddEntry(hnew, "New", "l");
605 leg->AddEntry(hold, "Old", "l");
606 leg->Draw();
607
608 c.cd(2);
609 gPad->SetGridy(1);
610 hratio->SetLineColor(kBlue);
611 hratio->SetStats(0);
612 hratio->SetTitle(Form("Ratio: new/old, %s;cos(#theta); New / Old", m_label[il].data()));
613 hratio->GetYaxis()->SetRangeUser(0.8, 1.2);
614 hratio->Draw("hist");
615
616 TLine* line = new TLine(m_cosMin, 1.0, m_cosMax, 1.0);
617 line->SetLineStyle(2);
618 line->Draw();
619
620 c.SaveAs(Form("cdcdedx_coscorr_fconsts_%s_%s.pdf", m_label[il].data(), m_suffix.data()));
621 c.SaveAs(Form("cdcdedx_coscorr_fconsts_%s_%s.root", m_label[il].data(), m_suffix.data()));
622
623 // cleanup
624 delete hnew;
625 delete hold;
626 delete hratio;
627 delete line;
628 }
629}
630
631//------------------------------------
633{
634
635 TCanvas cstats("cstats", "cstats", 1000, 500);
636 cstats.SetBatch(kTRUE);
637 cstats.Divide(2, 1);
638
639 cstats.cd(1);
640 auto hestats = getObjectPtr<TH1I>("hestats");
641 if (hestats) {
642 hestats->SetName(Form("hestats_%s", m_suffix.data()));
643 hestats->SetStats(0);
644 hestats->DrawCopy("");
645 }
646
647 cstats.cd(2);
648 auto htstats = getObjectPtr<TH1I>("htstats");
649 if (htstats) {
650 htstats->SetName(Form("htstats_%s", m_suffix.data()));
651 htstats->SetStats(0);
652 htstats->DrawCopy("");
653 }
654 cstats.Print(Form("cdcdedx_coscorr_stats_%s.pdf", m_suffix.data()));
655}
FitValues fitHistogram(TH1D *&hist)
Fit histogram with Gaussian and return mean, error, and width.
void getExpRunInfo()
function to extract calibration run/exp
TH1D * defineCosthHist(const std::string &tag)
function to define cosine histograms
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 plotCosThetaDist(TH1D *hCosth_all, TH1D *hCosth_pos, TH1D *hCosth_neg)
Plot cos(theta) distributions for all, positive, and negative tracks.
bool isMergePayload
merge payload at the time of calibration
unsigned int getRepresentativeLayer(unsigned int igroup) const
Representative CDC layer for each SL group (used to access group-wise constants): SL0 => 1,...
std::string m_suffix
add suffix to all plot name
DBObjPtr< CDCDedxCosineCor > m_DBCosineCor
Electron saturation correction DB object.
void plotConstants()
function to draw the old/new final constants
CDCDedxCosineAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
void fitGaussianWithRange(TH1D *&temphist, TString &status)
function to fit histogram in each cosine bin
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.
double m_sigLim
gaussian fit sigma limit
void plotFitResults(const std::vector< std::vector< double > > &dedxAll, const std::vector< std::vector< double > > &dedxNeg, const std::vector< std::vector< double > > &dedxPos)
Plot dE/dx fit results for all, negative, and positive tracks.
void plotdedxHist(std::vector< TH1D * > &hDedxCos_all, std::vector< TH1D * > &hDedxCos_neg, std::vector< TH1D * > &hDedxCos_pos)
function to draw the dE/dx histogram in costh bins
bool isMethodSep
if e+ e- need to be consider sep
void createPayload(std::vector< double > cosine)
function to store new payload after full calibration
double m_dedxMax
max dedx range for gain cal
bool isMakePlots
produce plots for status
double m_dedxMin
min dedx range for gain cal
void setHist(TH1D *h, int color, const char *title, double ymin, double ymax, int marker=20)
Set basic style (color, marker, title, y-range) for a TH1D histogram.
unsigned int m_cosBin
number of bins across cosine range
void defineHisto(std::vector< TH1D * > &hdedx, const std::string &tag, const std::string &chargeLabel)
function to define dE/dx histograms
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.
Container for Gaussian fit results of a histogram.
double meanErr
meanErr : uncertainty on the mean
double sigmaErr
sigmaErr : uncertainty on the width
TString status
status : fit status flag (e.g.