9#include <cdc/calibration/CDCdEdx/CDCDedxValidationAlgorithm.h>
10#include <cdc/calibration/CDCdEdx/CDCDedxWireGainAlgorithm.h>
12#include <cdc/dbobjects/CDCDedxWireGain.h>
13#include <cdc/dbobjects/CDCDedxCosineCor.h>
14#include <cdc/dbobjects/CDCDedx1DCell.h>
15#include <cdc/dbobjects/CDCDedxRunGain.h>
16#include <cdc/dbobjects/CDCDedxBadWires.h>
17#include <framework/database/IntervalOfValidity.h>
19#include <framework/database/Database.h>
20#include <framework/database/DBStore.h>
21#include <framework/database/Configuration.h>
50 m_eaMin(-TMath::Pi() / 2),
51 m_eaMax(+TMath::Pi() / 2),
67 std::vector<std::string> subdirs = {
"run",
"costh",
"mom",
"wire",
"injection",
"oneD"};
68 for (
const auto& dir : subdirs) {
69 gSystem->Exec(Form(
"mkdir -p plots/%s", dir.c_str()));
73 auto tBhabha = getObjectPtr<TTree>(
"tBhabha");
79 auto tRadee = getObjectPtr<TTree>(
"tRadee");
96 if (cruns == 0) B2INFO(
"start exp " << expRun.first <<
" and run " << expRun.second <<
"");
101 int estart = erStart.first;
102 int rstart = erStart.second;
107 else m_suffix = Form(
"e%d_r%d", estart, rstart);
113 auto ttree = getObjectPtr<TTree>(
"tRadee");
115 double dedx, costh, p, injtime = 0.0, injring = 1.0;
118 std::vector<double>* dedxhit = 0, *enta = 0;
119 std::vector<int>* layer = 0;
121 ttree->SetBranchAddress(
"dedx", &dedx);
122 ttree->SetBranchAddress(
"p", &p);
123 ttree->SetBranchAddress(
"costh", &costh);
124 ttree->SetBranchAddress(
"charge", &charge);
125 ttree->SetBranchAddress(
"injtime", &injtime);
126 ttree->SetBranchAddress(
"injring", &injring);
127 ttree->SetBranchAddress(
"dedxhit", &dedxhit);
128 ttree->SetBranchAddress(
"entaRS", &enta);
129 ttree->SetBranchAddress(
"layer", &layer);
131 std::vector<double> vtlocaledges;
133 m_tbins = vtlocaledges.size() - 1;
136 std::array<std::array<std::vector<TH1D*>, 2>, 13> hdedx_mom;
137 std::array<std::vector<TH1D*>, 2> hdedx_mom_peaks, hdedx_inj, hdedx_oned;
141 const double momBinW = (4.0 -
m_momMin) / 4;
143 std::string scos[13] = {
"acos",
"cos#theta > 0.0",
"cos#theta < 0.0",
"cos#theta <= -0.8",
144 "cos#theta > -0.8 and cos#theta <= -0.6",
145 "cos#theta > -0.6 and cos#theta <= -0.4",
"cos#theta > -0.4 and cos#theta <= -0.2",
146 "cos#theta > -0.2 and cos#theta <= 0",
"cos#theta > 0 and cos#theta <= 0.2",
147 "cos#theta > 0.2 and cos#theta <= 0.4",
"cos#theta > 0.4 and cos#theta <= 0.6",
148 "cos#theta > 0.6 and cos#theta <= 0.8",
"cos#theta > 0.8"
150 std::string stype[2] = {
"posi",
"elec"};
151 std::string sLayer[2] = {
"IL",
"OL"};
154 for (
int ic = 0; ic < 13; ic++) {
155 for (
int it = 0; it < 2; ++it) {
157 defineHisto(hdedx_mom[ic][it],
"mom", Form(
"%d_%s", ic, stype[it].data()));
162 for (
unsigned int ir = 0; ir < 2; ir++) {
164 hdedx_mom_peaks[ir].resize(4);
165 hdedx_oned[ir].resize(
m_eaBin);
168 defineHisto(hdedx_mom_peaks[ir],
"mom_peaks", Form(
"%s", stype[ir].data()));
169 defineHisto(hdedx_oned[ir],
"oned", Form(
"%s", sLayer[ir].data()));
173 double icos[3] = {0, -1, -1};
177 for (
int i = 0; i < ttree->GetEntries(); ++i) {
181 if (dedx <= 0 || injtime < 0 || injring < 0)
continue;
184 int binIndex =
static_cast<int>((abs(p) -
m_momMin) / momBinWidth);
188 icos[1] = (costh > 0) ? 1 : 2;
189 icos[2] = int((costh + 1.0) / 0.2) + 3;
190 if (icos[2] < 3) icos[2] = 3;
191 if (icos[2] > 12) icos[2] = 12;
194 chgtype = (charge > 0) ? 0 : 1;
197 if (binIndex >= 0 && binIndex <
m_momBins) {
198 hdedx_mom[icos[0]][chgtype][binIndex]->Fill(dedx);
199 hdedx_mom[icos[1]][chgtype][binIndex]->Fill(dedx);
200 hdedx_mom[icos[2]][chgtype][binIndex]->Fill(dedx);
207 int wr = (injring > 0.5) ? 1 : 0;
210 unsigned int tb = htimes->GetXaxis()->FindBin(injtime);
211 tb = std::min(tb,
static_cast<unsigned int>(
m_tbins)) - 1;
214 htimes->Fill(injtime);
215 hdedx_inj[wr][tb]->Fill(dedx);
218 int binI =
static_cast<int>((abs(p) -
m_momMin) / momBinW);
219 if (binI >= 0 && binI < 4) {
220 hdedx_mom_peaks[chgtype][binI]->Fill(dedx);
224 for (
unsigned int j = 0; j < dedxhit->size(); ++j) {
225 if (dedxhit->at(j) == 0)
continue;
227 double entaval = enta->at(j);
228 int ibin = std::floor((entaval -
m_eaMin) / eaBW);
229 if (ibin < 0 || ibin >=
m_eaBin)
continue;
231 int mL = (layer->at(j) < 8) ? 0 : 1;
232 hdedx_oned[mL][ibin]->Fill(dedxhit->at(j));
237 for (
int ic = 0; ic < 13; ic++) {
238 for (
int it = 0; it < 2; ++it) {
239 printCanvas(hdedx_mom[ic][it], Form(
"plots/mom/dedx_vs_mom_%d_%s_%s", ic, stype[it].data(),
m_suffix.data()),
"mom");
242 for (
int it = 0; it < 2; ++it) {
244 printCanvas(hdedx_oned[it], Form(
"plots/oneD/dedx_vs_1D_%s_%s", sLayer[it].data(),
m_suffix.data()),
"oned");
253 auto ttree = getObjectPtr<TTree>(
"tBhabha");
258 std::vector<int>* wire = 0;
259 ttree->SetBranchAddress(
"wire", &wire);
261 std::vector<double>* dedxhit = 0;
262 ttree->SetBranchAddress(
"dedxhit", &dedxhit);
264 ttree->SetBranchAddress(
"dedx", &dedx);
265 ttree->SetBranchAddress(
"run", &run);
266 ttree->SetBranchAddress(
"charge", &charge);
267 ttree->SetBranchAddress(
"costh", &costh);
269 std::map<int, TH1D*> hdedx_run;
270 std::array<std::vector<TH1D*>, 3> hdedx_cos;
271 std::array<std::vector<TH1D*>, 2> hdedx_cos_peaks;
272 std::vector<TH1D*> hdedxhit(c_nSenseWires);
277 std::string stype[3] = {
"all",
"posi",
"elec"};
279 for (
int it = 0; it < 3; ++it) {
284 for (
int ir = 0; ir < 2; ir++) {
285 hdedx_cos_peaks[ir].resize(4);
286 defineHisto(hdedx_cos_peaks[ir],
"cos_peaks", Form(
"%s", stype[ir + 1].data()));
292 for (
int i = 0; i < ttree->GetEntries(); ++i) {
294 if (dedx <= 0)
continue;
297 if (hdedx_run.find(run) == hdedx_run.end()) {
298 std::string histName = Form(
"hist_dedx_run_%d", run);
299 std::string histTitle = Form(
"dE/dx Histogram for Run %d", run);
304 hdedx_run[run]->Fill(dedx);
307 int binIndex =
static_cast<int>((costh -
m_cosMin) / cosBinWidth);
308 if (binIndex >= 0 && binIndex <
m_cosBins) {
309 hdedx_cos[0][binIndex]->Fill(dedx);
312 hdedx_cos[1][binIndex]->Fill(dedx);
314 hdedx_cos[2][binIndex]->Fill(dedx);
318 for (
unsigned int j = 0; j < wire->size(); ++j) {
319 int jwire = wire->at(j);
320 double jhitdedx = dedxhit->at(j);
321 hdedxhit[jwire]->Fill(jhitdedx);
325 int binI =
static_cast<int>((costh -
m_cosMin) / cosBinW);
326 if (binI >= 0 && binI < 4) {
328 hdedx_cos_peaks[0][binI]->Fill(dedx);
330 hdedx_cos_peaks[1][binI]->Fill(dedx);
335 printCanvas(hdedx_cos[0], Form(
"plots/costh/dedx_vs_cos_all_%s",
m_suffix.data()),
"costh");
336 printCanvas(hdedx_cos[1], Form(
"plots/costh/dedx_vs_cos_positrons_%s",
m_suffix.data()),
"costh");
337 printCanvas(hdedx_cos[2], Form(
"plots/costh/dedx_vs_cos_electrons_%s",
m_suffix.data()),
"costh");
347 double binWidth = 0.0;
351 }
else if (var ==
"oned") {
353 }
else if (var ==
"costh") {
355 }
else if (var ==
"inj") {
357 }
else if (var ==
"mom_peaks") {
358 xbins = 4; xmin =
m_momMin; xmax = 4.0;
359 }
else if (var ==
"cos_peaks") {
365 if (var ==
"costh" || var ==
"mom" || var ==
"mom_peaks" || var ==
"cos_peaks" || var ==
"oned") {
366 binWidth = (xmax - xmin) / xbins;
369 for (
int ic = 0; ic < xbins; ic++) {
370 std::string title = Form(
"dedxhit-dist, wire:%d", ic);
371 std::string name = Form(
"hdedx_%s_%d", var.data(), ic);
373 if (var ==
"costh" || var ==
"mom" || var ==
"mom_peaks" || var ==
"cos_peaks" || var ==
"oned") {
374 double min = ic * binWidth + xmin;
375 double max = min + binWidth;
376 title = Form(
"%s: (%0.02f, %0.02f) %s", var.data(), min, max, stype.data());
377 name = Form(
"hdedx_%s_%s_%d", var.data(), stype.data(), ic);
378 }
else if (var ==
"inj") {
380 title = Form(
"%s, time(%s)", stype.data(), label.data());
381 name = Form(
"h%s_%s_t%d", var.data(), stype.data(), ic);
384 htemp[ic]->SetTitle(Form(
"%s;dedx;entries", title.data()));
395 }
else if (svar ==
"costh") {
398 B2FATAL(
"wrong input");
400 double binWidth = (xmax - xmin) / xbins;
403 TCanvas* ctmp =
new TCanvas(
"tmp",
"tmp", 1200, 1200);
407 std::stringstream psname;
408 psname << Form(
"%s.pdf[", namesfx.data());
409 ctmp->Print(psname.str().c_str());
411 psname << Form(
"%s.pdf", namesfx.data());
414 for (
int i = 0; i < xbins; ++i) {
418 double emean, emeanErr, esigma, esigmaErr;
419 double pmean, pmeanErr, psigma, psigmaErr;
421 fit(htemp[0][i], emean, emeanErr, esigma, esigmaErr);
422 fit(htemp[1][i], pmean, pmeanErr, psigma, psigmaErr);
424 double min = i * binWidth + xmin;
425 double max = min + binWidth;
427 TPaveText pt(0.6, 0.63, 0.85, 0.89,
"NBNDC");
430 pt.AddText(Form(
"#mu_{fit}: %0.03f#pm%0.03f", emean, emeanErr));
431 pt.AddText(Form(
"#sigma_{fit}: %0.03f#pm%0.03f", esigma, esigmaErr));
434 pt.AddText(Form(
"#mu_{fit}: %0.03f#pm%0.03f", pmean, pmeanErr));
435 pt.AddText(Form(
"#sigma_{fit}: %0.03f#pm%0.03f", psigma, psigmaErr));
437 htemp[0][i]->SetStats(0);
438 htemp[1][i]->SetStats(0);
439 htemp[0][i]->SetFillColor(0);
440 htemp[1][i]->SetFillColor(0);
441 htemp[0][i]->SetLineColor(8);
442 htemp[1][i]->SetLineColor(9);
443 htemp[0][i]->SetTitle(Form(
"%s: (%0.02f, %0.02f)", svar.data(), min, max));
444 if (htemp[0][i]->GetEntries() > 0)
445 htemp[0][i]->Scale(1.0 / htemp[0][i]->GetEntries());
446 if (htemp[1][i]->GetEntries() > 0)
447 htemp[1][i]->Scale(1.0 / htemp[1][i]->GetEntries());
449 if (htemp[1][i]->GetMaximum() > htemp[0][i]->GetMaximum())
450 htemp[0][i]->SetMaximum(htemp[1][i]->GetMaximum());
452 htemp[0][i]->DrawCopy(
"HIST");
453 htemp[1][i]->DrawCopy(
"same HIST");
454 pt.DrawClone(
"same");
456 TLegend* lego =
new TLegend(0.15, 0.67, 0.3, 0.8);
457 lego->AddEntry(htemp[0][i],
"e+",
"l");
458 lego->AddEntry(htemp[1][i],
"e-",
"l");
461 if ((i + 1) % 4 == 0 || i == xbins - 1) {
462 ctmp->SetBatch(kTRUE);
463 ctmp->Print(psname.str().c_str());
464 if ((i + 1) % 4 == 0) ctmp->Clear(
"D");
469 psname << Form(
"%s.pdf]", namesfx.data());
470 ctmp->Print(psname.str().c_str());
482 }
else if (svar ==
"oned") {
484 }
else if (svar ==
"costh") {
486 }
else if (svar ==
"inj") {
488 }
else if (svar ==
"mom_peaks") {
489 xbins = 4; xmin =
m_momMin; xmax = 4.0;
491 B2FATAL(
"wrong input");
495 TCanvas* ctmp =
new TCanvas(
"tmp",
"tmp", 1200, 1200);
499 std::stringstream psname;
500 psname << Form(
"%s.pdf[", namesfx.data());
501 ctmp->Print(psname.str().c_str());
503 psname << Form(
"%s.pdf", namesfx.data());
505 std::ofstream outFile;
506 outFile.open(Form(
"%s.txt", namesfx.data()));
510 for (
int i = 0; i < xbins; ++i) {
512 ctmp->cd(i % 16 + 1);
513 TPaveText pt(0.6, 0.73, 0.85, 0.89,
"NBNDC");
516 if (svar ==
"oned") {
517 unsigned int minbin, maxbin;
519 htemp[i]->SetTitle(Form(
"dedxhit-dist, entabin: %d ;%d;%d", i, minbin, maxbin));
523 const double binWidth = (xmax - xmin) / xbins;
524 double binCenter = xmin + (i + 0.5) * binWidth;
526 outFile << binCenter <<
" " << dedxmean << std::endl;
528 double mean, meanErr, sigma, sigmaErr;
529 fit(htemp[i], mean, meanErr, sigma, sigmaErr);
531 if (svar ==
"mom" || svar ==
"costh" || svar ==
"mom_peaks") {
532 const double binWidth = (xmax - xmin) / xbins;
533 double binCenter = xmin + (i + 0.5) * binWidth;
535 outFile << binCenter <<
" " << mean <<
" " << meanErr <<
" " << sigma <<
" " << sigmaErr << std::endl;
538 outFile << i <<
" " << label <<
" " << mean <<
" " << meanErr <<
" " << sigma <<
" " << sigmaErr << std::endl;
541 pt.AddText(Form(
"#mu_{fit}: %0.03f#pm%0.03f", mean, meanErr));
542 pt.AddText(Form(
"#sigma_{fit}: %0.03f#pm%0.03f", sigma, sigmaErr));
544 htemp[i]->SetStats(0);
545 htemp[i]->DrawCopy(
"");
546 pt.DrawClone(
"same");
548 if ((i + 1) % 16 == 0 || ((i + 1) == xbins)) {
549 ctmp->SetBatch(kTRUE);
550 ctmp->Print(psname.str().c_str());
557 psname << Form(
"%s.pdf]", namesfx.data());
558 ctmp->Print(psname.str().c_str());
568 std::string status =
"";
570 if (hist->Integral() > 100)
573 if (status !=
"fitOK") {
574 hist->SetFillColor(kOrange);
575 mean = 0.0, meanErr = 0.0, sigma = 0.0, sigmaErr = 0.0;
577 mean = hist->GetFunction(
"gaus")->GetParameter(1);
578 meanErr = hist->GetFunction(
"gaus")->GetParError(1);
579 sigma = hist->GetFunction(
"gaus")->GetParameter(2);
580 sigmaErr = hist->GetFunction(
"gaus")->GetParError(2);
581 hist->SetFillColor(kYellow);
588 TCanvas* ctmp =
new TCanvas(
"tmp",
"tmp", 1200, 1200);
592 std::stringstream psname;
593 psname << Form(
"%s.pdf[", namesfx.data());
594 ctmp->Print(psname.str().c_str());
596 psname << Form(
"%s.pdf", namesfx.data());
598 std::ofstream outFile;
599 outFile.open(Form(
"%s.txt", namesfx.data()));
603 for (
const auto& entry : htemp) {
604 int run = entry.first;
605 TH1D* hist = entry.second;
607 ctmp->cd(irun % 16 + 1);
609 TPaveText pt(0.6, 0.73, 0.85, 0.89,
"NBNDC");
612 double mean, meanErr, sigma, sigmaErr;
613 fit(hist, mean, meanErr, sigma, sigmaErr);
615 outFile << run <<
" " << mean <<
" " << meanErr <<
" " << sigma <<
" " << sigmaErr << std::endl;
617 pt.AddText(Form(
"#mu_{fit}: %0.03f#pm%0.03f", mean, meanErr));
618 pt.AddText(Form(
"#sigma_{fit}: %0.03f#pm%0.03f", sigma, sigmaErr));
622 pt.DrawClone(
"same");
624 if ((irun + 1) % 16 == 0 || irun ==
int(htemp.size() - 1)) {
625 ctmp->SetBatch(kTRUE);
626 ctmp->Print(psname.str().c_str());
632 ctmp->Print(psname.str().c_str());
634 psname << Form(
"%s.pdf]", namesfx.data());
635 ctmp->Print(psname.str().c_str());
645 double histmean = temphist->GetMean();
646 double histrms = temphist->GetRMS();
647 temphist->GetXaxis()->SetRangeUser(histmean - 5.0 * histrms, histmean + 5.0 * histrms);
649 int fs = temphist->Fit(
"gaus",
"Q0");
651 B2INFO(Form(
"\tFit (round 1) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
652 status =
"fitFailed";
655 double mean = temphist->GetFunction(
"gaus")->GetParameter(1);
656 double width = temphist->GetFunction(
"gaus")->GetParameter(2);
657 temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
658 fs = temphist->Fit(
"gaus",
"QR",
"", mean -
m_sigmaR * width, mean +
m_sigmaR * width);
660 B2INFO(Form(
"\tFit (round 2) for hist (%s) failed (status = %d)", temphist->GetName(), fs));
661 status =
"fitFailed";
664 temphist->GetXaxis()->SetRangeUser(mean - 5.0 * width, mean + 5.0 * width);
665 B2INFO(Form(
"\tFit for hist (%s) successful (status = %d)", temphist->GetName(), fs));
674 std::vector<double> vdedx_means;
675 std::vector<double> layermean(c_maxNSenseLayers);
676 std::vector<double> lgmean(c_maxNSenseLayers);
678 std::ofstream outFile, outFileLayer, outFileAvg, outFilebdwire;
679 outFile.open(Form(
"plots/wire/dedx_mean_gwire_%s.txt",
m_suffix.data()));
680 outFilebdwire.open(Form(
"plots/wire/dedx_mean_badwire_%s.txt",
m_suffix.data()));
681 outFileLayer.open(Form(
"plots/wire/dedx_mean_layer_%s.txt",
m_suffix.data()));
682 outFileAvg.open(Form(
"plots/wire/dedx_mean_layer_avg_%s.txt",
m_suffix.data()));
684 int activelayers = 0;
685 double layeravg = 0.0;
693 for (
unsigned int il = 0; il < c_maxNSenseLayers; ++il) {
695 int activewires = 0, goodwires = 0;
699 for (
unsigned int iw = 0; iw < cdcgeo.
nWiresInLayer(il); ++iw) {
702 unsigned int minbin, maxbin;
704 hdedxhit[jwire]->SetTitle(Form(
"dedxhit-dist, wire: %d ;%d;%d", jwire, minbin, maxbin));
707 vdedx_means.push_back(dedxmean);
709 if (Badwire->getBadWireStatus(jwire) == kTRUE)
710 outFilebdwire << jwire <<
" " << dedxmean << std::endl;
712 outFile << jwire <<
" " << dedxmean << std::endl;
714 if (vdedx_means.at(jwire) > 0) {
715 layermean[il] += vdedx_means.at(jwire);
717 if (Badwire->getBadWireStatus(jwire) != kTRUE) {
718 lgmean[il] += vdedx_means.at(jwire);
724 if (activewires > 0) layermean[il] /= activewires;
725 else layermean[il] = 1.0;
727 if (goodwires > 0) lgmean[il] /= goodwires;
728 else lgmean[il] = 1.0;
730 outFileLayer << il <<
" " << layermean[il] <<
" " << lgmean[il] << std::endl;
733 if (il >= 8 && layermean[il] > 0) {
734 layeravg += layermean[il];
740 if (activelayers > 0) layeravg /= activelayers;
741 outFileAvg << layeravg << std::endl;
744 outFilebdwire.close();
745 outFileLayer.close();
751 const std::vector<double>& vdedx_mean)
753 TCanvas* ctmp =
new TCanvas(
"tmp",
"tmp", 900, 900);
756 std::stringstream psname;
757 psname << Form(
"%s.pdf[", namesfx.data());
758 ctmp->Print(psname.str().c_str());
760 psname << Form(
"%s.pdf", namesfx.data());
762 for (
unsigned int ip = 0; ip < c_nwireCDC; ip++) {
763 int minbin = std::stoi(temp[ip]->GetXaxis()->GetTitle());
764 int maxbin = std::stoi(temp[ip]->GetYaxis()->GetTitle());
765 temp[ip]->SetFillColor(kYellow - 9);
766 temp[ip]->SetTitle(Form(
"%s, #mu_{trunc} %0.03f;dedxhit;entries", temp[ip]->GetTitle(), vdedx_mean.at(ip)));
768 ctmp->cd(ip % 16 + 1);
770 temp[ip]->DrawCopy(
"hist");
771 TH1D* hdedxhitC = (TH1D*)temp[ip]->Clone(Form(
"%sC", temp[ip]->GetName()));
772 hdedxhitC->GetXaxis()->SetRange(minbin, maxbin);
773 hdedxhitC->SetFillColor(kAzure + 1);
774 hdedxhitC->DrawCopy(
"same histo");
776 if ((ip + 1) % 16 == 0) {
777 ctmp->SetBatch(kTRUE);
778 ctmp->Print(psname.str().c_str());
787 psname << Form(
"%s.pdf]", namesfx.data());
788 ctmp->Print(psname.str().c_str());
795 for (
int ib = 0; ib < 69; ib++) {
796 fixedges[ib] = ib * 0.5 * 1e3;
797 if (ib > 40 && ib <= 60) fixedges[ib] = fixedges[ib - 1] + 1.0 * 1e3;
798 else if (ib > 60 && ib <= 64) fixedges[ib] = fixedges[ib - 1] + 10.0 * 1e3;
799 else if (ib > 64 && ib <= 65) fixedges[ib] = fixedges[ib - 1] + 420.0 * 1e3;
800 else if (ib > 65 && ib <= 66) fixedges[ib] = fixedges[ib - 1] + 500.0 * 1e3;
801 else if (ib > 66) fixedges[ib] = fixedges[ib - 1] + 2e6;
802 vtlocaledges.push_back(fixedges[ib]);
810 TCanvas cstats(
"cstats",
"cstats", 800, 400);
811 cstats.SetBatch(kTRUE);
815 auto hestats = getObjectPtr<TH1I>(
"hestats");
817 hestats->SetName(Form(
"hestats_%s",
m_suffix.data()));
818 hestats->SetStats(0);
819 hestats->DrawCopy(
"");
823 auto htstats = getObjectPtr<TH1I>(
"htstats");
825 htstats->SetName(Form(
"htstats_%s",
m_suffix.data()));
826 htstats->SetStats(0);
827 htstats->DrawCopy(
"");
830 cstats.Print(Form(
"cdcdedx_stats_%s.pdf",
m_suffix.data()));
841 dbConfiguration.setGlobalTags({
"online"});
847 B2FATAL(
"Setting both testing payload and Global Tag or setting no one of them.");
867 std::vector<double> wiregain;
868 std::vector<double> layermean(c_maxNSenseLayers);
871 if (!DBWireGains.
isValid()) B2FATAL(
"Wire gain data are not valid.");
876 for (
unsigned int il = 0; il < c_maxNSenseLayers; ++il) {
881 for (
unsigned int iw = 0; iw < cdcgeo.
nWiresInLayer(il); ++iw) {
884 wiregain.push_back(DBWireGains->getWireGain(jwire));
886 if (wiregain.at(jwire) > 0) {
887 layermean[il] += wiregain.at(jwire);
892 if (activewires > 0) layermean[il] /= activewires;
893 else layermean[il] = 1.0;
897 return { wiregain, layermean };
905 std::vector<double> cosgain, cos;
908 if (!DBCosineCor.
isValid()) B2FATAL(
"Cosine gain data are not valid.");
910 unsigned int nCosBins = DBCosineCor->getSize();
912 for (
unsigned int il = 0; il < nCosBins; ++il) {
913 double costh = -1.0 + (il + 0.5) * 2.0 / nCosBins;
915 cosgain.push_back(DBCosineCor->getMean(il));
916 cos.push_back(costh);
920 return {cosgain, cos};
928 std::vector<double> inner1D, outer1D, Enta;
931 if (!DBOneDCell.
isValid()) B2FATAL(
"OneD cell gain data are not valid.");
933 for (
int i = 0; i < 2; i++) {
935 unsigned int nBins = DBOneDCell->getNBins(i);
936 double binSize = TMath::Pi() / nBins;
938 for (
unsigned int nbin = 0; nbin < nBins; nbin++) {
940 double enta = (-1 * TMath::Pi() / 2.0) + binSize * nbin;
942 Enta.push_back(enta);
943 inner1D.push_back(DBOneDCell->getMean(0, nbin));
945 outer1D.push_back(DBOneDCell->getMean(17, nbin));
949 return {inner1D, outer1D, Enta};
959 if (!RunGain.
isValid()) B2FATAL(
"Run gain data are not valid.");
960 double gain = RunGain->getRunGain();
void bhabhaValidation()
Validate dE/dx using bhabha sample (vs run, cosine)
double m_eaMax
upper edge of entrance angle
double m_momMin
min range of momentum
void resetDatabase()
Clear current DB pointers and state.
void wireGain(std::vector< TH1D * > &hdedxhit)
Validate wire gain data using dE/dx histograms.
void printCanvasdEdx(std::array< std::vector< TH1D * >, 2 > &htemp, std::string namesfx, std::string svar)
Draw dE/dx histograms for momentum and cosine bins.
void printCanvasRun(std::map< int, TH1D * > &htemp, std::string namesfx)
Draw dE/dx per run histogram canvas.
int m_cosBins
bins for cosine
void radeeValidation()
Validate dE/dx using radee sample (vs momentum, injection time)
double m_sigmaR
fit dedx dist in sigma range
void setTextCosmetics(TPaveText pt, Color_t color)
Set text cosmetics for TPaveText.
WireGainData getwiregain(int experiment, int run)
Retrieve wire gain data from DB.
void getExpRunInfo()
function to get extract calibration run/exp
double m_cosMax
max range of cosine
OnedData getonedgain(int experiment, int run)
Retrieve 1D gain data from DB.
double * m_tedges
internal time array (copy of vtlocaledges)
std::string m_GlobalTagName
Global Tag name.
void DatabaseIN(int experiment, int run)
Load database payload for given run.
void defineTimeBins(std::vector< double > &vtlocaledges)
Set bin edges for injection time.
std::array< std::string, 2 > m_sring
injection ring name
CosGainData getcosgain(int experiment, int run)
Retrieve cosine gain data from DB.
void printCanvasWire(std::vector< TH1D * > temp, std::string namesfx, const std::vector< double > &vdedx_mean)
Plot dE/dx vs wire number.
std::string m_testingPayloadName
Testing payload location.
std::string m_suffix
suffix string to separate plots
int m_momBins
bins for momentum
int m_dedxBins
bins for dedx histogram
DBObjPtr< CDCGeometry > m_cdcGeo
Geometry of CDC.
double m_momMax
max range of momentum
CDCDedxValidationAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
double getrungain(int experiment, int run)
Retrieve run gain data from DB.
void fitGaussianWRange(TH1D *&temphist, std::string &status)
Perform Gaussian fit with range on a histogram.
double m_cosMin
min range of cosine
void plotEventStats()
Plot summary statistics of selected events.
virtual EResult calibrate() override
Main calibration method.
void printCanvas(std::vector< TH1D * > &htemp, std::string namesfx, std::string svar)
Draw dE/dx histograms across bins.
double m_dedxMax
max range of dedx
StoreObjPtr< EventMetaData > m_EventMetaData
Event metadata.
double m_dedxMin
min range of dedx
unsigned int m_tbins
internal time bins
void defineHisto(std::vector< TH1D * > &htemp, std::string var, std::string stype)
Define dE/dx histograms for plotting.
std::string getTimeBinLabel(const double &tedges, const int &it)
Get time bin label string.
void fit(TH1D *&hist, double &mean, double &meanErr, double &sigma, double &sigmaErr)
Perform full Gaussian fit and extract parameters.
double m_eaMin
lower edge of entrance angle
A calibration algorithm for CDC dE/dx wire gains.
void getTruncatedBins(TH1D *hdedxhit, unsigned int &binlow, unsigned int &binhigh)
function to get bins of trunction from histogram
double getTruncationMean(TH1D *hdedxhit, int binlow, int binhigh)
function to get mean of trunction from histogram
The Class for CDC Geometry Parameters.
unsigned nWiresInLayer(int layerId) const
Returns wire numbers in a layer.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
Base class for calibration algorithms.
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 successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
static Configuration & getInstance()
Get a reference to the instance which will be used when the Database is initialized.
bool isValid() const
Check whether a valid object was obtained from the database.
Class for accessing objects in the database.
Singleton class to cache database objects.
static DataStore & Instance()
Instance of singleton Store.
void setInitializeActive(bool active)
Setter for m_initializeActive.
void reset(EDurability durability)
Frees memory occupied by data store items and removes all objects from the map.
void reset(bool keepEntries=false)
Invalidate all payloads.
static Database & Instance()
Instance of a singleton Database.
static DBStore & Instance()
Instance of a singleton DBStore.
void updateEvent()
Updates all intra-run dependent objects.
void update()
Updates all objects that are outside their interval of validity.
static void reset(bool keepConfig=false)
Reset the database instance.
Abstract base class for different kinds of events.
Container for cosine gain data.
Container for 1D gain data.
Container for wire gain data.