9#include <cdc/calibration/CDCdEdx/CDCDedx1DCellAlgorithm.h>
19 m_eaMin(-TMath::Pi() / 2),
20 m_eaMax(+TMath::Pi() / 2),
41 setDescription(
"A calibration algorithm for the CDC dE/dx entrance angle cleanup correction");
53 B2FATAL(
"There is no valid previous payload for CDCDedx1DCell");
56 auto ttree = getObjectPtr<TTree>(
"tree");
59 std::vector<double>* dedxhit = 0, *enta = 0;
60 std::vector<int>* layer = 0;
61 double pt = 0, costh = 0;
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);
82 std::vector<TH1D*> hdedxhit[2];
86 TH2D* hptcosth =
new TH2D(
"hptcosth",
"pt vs costh dist;pt;costh", 1000, -8.0, 8.0, 1000, -1.0, 1.0);
91 for (
int i = 0; i < ttree->GetEntries(); ++i) {
95 if (std::abs(costh) >
m_cosMax)
continue;
101 if (std::abs(pt) >
m_ptMax)
continue;
104 int rand = gRandom->Integer(100);
105 if (rand < 10) hptcosth->Fill(pt, costh);
107 for (
unsigned int j = 0; j < dedxhit->size(); ++j) {
109 if (dedxhit->at(j) == 0)
continue;
111 double entaval = enta->at(j);
114 if (ibin < 0 || ibin >
m_eaBin)
continue;
117 if (layer->at(j) < 8)mL = 0;
120 hdedxlay[mL]->Fill(dedxhit->at(j));
121 if (rand < 10) hentalay[mL]->Fill(entaval);
125 hdedxhit[mL][jbinea]->Fill(dedxhit->at(j));
129 for (
int il = 0; il < 2; il++) {
131 int minlay = 0, maxlay = 0;
135 hdedxlay[il]->SetTitle(Form(
"%s;%d;%d", hdedxlay[il]->GetTitle(), minlay, maxlay));
138 std::vector<double>tempconst;
148 TH1D* htemp = (TH1D*)hdedxhit[il][jea]->Clone(Form(
"h_%s_b%d_c",
m_label[il].data(), jea));
150 int minbin = 1, maxbin = 1;
161 tempconst.push_back(dedxmean);
163 hdedxhit[il][iea]->SetTitle(Form(
"%s, #mu_{truc} = %0.5f;%d;%d", hdedxhit[il][iea]->GetTitle(), dedxmean, minbin, maxbin));
167 std::vector<double>layerconst;
170 for (
int iea = 0; iea <
m_eaBin; iea++) {
173 layerconst.push_back(tempconst.at(jea));
205 for (
int il = 0; il < 2; il++) {
209 delete hdedxhit[il][iea];
214 for (
int il = 0; il < 2; il++) {
227 if (cruns == 0) B2INFO(
"CDCDedxBadWires: start exp " << expRun.first <<
" and run " << expRun.second <<
"");
232 int estart = erStart.first;
233 int rstart = erStart.second;
236 int eend = erEnd.first;
237 int rend = erEnd.second;
241 m_runExp = Form(
"Range (%d:%d,%d:%d)", estart, rstart, eend, rend);
243 else m_suffix = Form(
"e%d_r%dr%d", estart, rstart, rend);
250 std::map<int, std::vector<double>> bounds;
251 std::map<int, std::vector<int>> steps;
253 const std::array<int, 2> nDev{8, 4};
254 bounds[0] = {0, 108, 123, 133, 158, 183, 193, 208, 316};
255 steps[0] = {9, 3, 2, 1, 1, 2, 3, 9};
256 bounds[1] = {0, 38, 158, 278, 316};
257 steps[1] = {2, 1, 1, 2};
259 for (
int il = 0; il < 2; il++) {
261 for (
int ibin = 0; ibin <= nDev[il]; ibin++) bounds[il][ibin] = bounds[il][ibin] *
m_binSplit;
263 int ieaprime = -1, temp = -99, ibin = 0;
268 for (
int iea = 0; iea <
m_eaBin; iea++) {
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;
278 if (ieaprime != temp) {
281 double binvalue = pastbin + binwidth;
283 if (std::abs(binvalue) < 1e-5)binvalue = 0;
296 for (
int il = 0; il < 2; il++) {
298 std::string title = Form(
"dedxhit dist (%s): %s ; dedxhit;entries",
m_label[il].data(),
m_runExp.data());
300 hdedxlay[il]->SetTitle(Form(
"%s", title.data()));
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());
308 hentalay[il] =
new TH1D(Form(
"hentalay%s",
m_label[il].data()),
"",
m_eaBinLocal[il], nvarBins);
309 hentalay[il]->SetTitle(Form(
"%s", title.data()));
315 double width = max - min;
317 if (
isPrintLog) B2INFO(
"bin: " << iea <<
" ], min:" << min <<
" , max: " << max <<
" , width: " << width);
319 title = Form(
"%s: entaRS = (%0.03f to %0.03f)",
m_label[il].data(), min, max);
321 hdedxhit[il][iea]->SetTitle(Form(
"%s", title.data()));
331 double sum = hist->Integral();
332 if (sum <= 0 || hist->GetNbinsX() <= 0) {
333 binlow = 1; binhigh = 1;
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);
358 if (hist->Integral() < 100)
return 1.0;
360 if (binlow <= 0 || binhigh > hist->GetNbinsX())
return 1.0;
362 double binweights = 0., sumofbc = 0.;
363 for (
int ibin = binlow; ibin <= binhigh; ibin++) {
364 double bcdedx = hist->GetBinContent(ibin);
366 binweights += (bcdedx * hist->GetBinCenter(ibin));
370 if (sumofbc > 0)
return binweights / sumofbc;
378 B2INFO(
"dE/dx one cell calibration: Generating payloads");
380 for (
unsigned int il = 0; il < 2; il++) {
385 B2ERROR(
"merging failed because of unmatch bins (old " <<
m_eaBin <<
" new " << nbins <<
")");
387 for (
unsigned int iea = 0; iea < nbins; iea++) {
401 B2INFO(
"dE/dx Calibration done for CDCDedx1DCell");
408 std::map<
int, std::vector<int>> steps)
411 TCanvas cmfactor(
"cmfactor",
"Merging factors", 800, 400);
412 cmfactor.Divide(2, 1);
414 for (
int il = 0; il < 2; il++) {
416 nvarBins = &bounds[il][0];
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()));
421 for (
int ibin = 0; ibin < nDev[il]; ibin++) hist->SetBinContent(ibin + 1, steps[il][ibin]);
424 hist->SetFillColor(kYellow);
429 cmfactor.SaveAs(Form(
"cdcdedx_1dcell_mergefactor%s.pdf",
m_suffix.data()));
430 cmfactor.SaveAs(Form(
"cdcdedx_1dcell_mergefactor%s.root",
m_suffix.data()));
437 TCanvas ctmp(
"tmp",
"tmp", 1200, 1200);
439 std::stringstream psname;
441 psname << Form(
"cdcdedx_1dcell_dedxhit%s.pdf[",
m_suffix.data());
442 ctmp.Print(psname.str().c_str());
444 psname << Form(
"cdcdedx_1dcell_dedxhit%s.pdf",
m_suffix.data());
446 for (
int il = 0; il < 2; il++) {
450 int minbin = std::stoi(hdedxhit[il][jea]->GetXaxis()->GetTitle());
451 int maxbin = std::stoi(hdedxhit[il][jea]->GetYaxis()->GetTitle());
453 ctmp.cd(jea % 16 + 1);
454 hdedxhit[il][jea]->SetFillColor(4 + il);
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");
464 ctmp.Print(psname.str().c_str());
472 psname << Form(
"cdcdedx_1dcell_dedxhit%s.pdf]",
m_suffix.data());
473 ctmp.Print(psname.str().c_str());
480 TCanvas cdedxlayer(
"layerdedxhit",
"Inner and Outer Layer dedxhit dist", 900, 400);
481 cdedxlayer.Divide(2, 1);
483 for (
int il = 0; il < 2; il++) {
484 int minlay = 0, maxlay = 0;
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));
493 cdedxlayer.cd(il + 1);
494 hdedxlay[il]->SetFillColor(kYellow);
495 hdedxlay[il]->Draw(
"histo");
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");
505 cdedxlayer.SaveAs(Form(
"cdcdedx_1dcell_dedxlay%s.pdf",
m_suffix.data()));
506 cdedxlayer.SaveAs(Form(
"cdcdedx_1dcell_dedxlay%s.root",
m_suffix.data()));
513 TCanvas ceadist(
"ceadist",
"Enta distributions", 800, 400);
514 ceadist.Divide(2, 1);
516 for (
int il = 0; il < 2; il++) {
520 hentalay[il]->SetFillColor(kYellow);
521 hentalay[il]->Draw(
"hist");
524 TCanvas cptcos(
"cptcos",
"pt vs costh dist.", 400, 400);
526 hptcosth->Draw(
"colz");
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()));
537 TH1D* hconst, *hconstvar;
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()));
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()));
552 hconstvar->SetBinContent(iea + 1, tempconst.at(iea));
555 for (
int jea = 0; jea <
m_eaBin; jea++) hconst->SetBinContent(jea + 1, layerconst.at(jea));
557 gStyle->SetOptStat(
"ne");
558 TCanvas cconst(
"cconst",
"Calirbation Constants", 800, 400);
561 cconst.SetWindowSize(1000, 800);
565 hconst->SetFillColor(kYellow);
566 hconst->Draw(
"histo");
569 hconstvar->SetFillColor(kBlue);
570 hconstvar->Draw(
"hist");
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()));
584 TH1D* hnewconst[2], *holdconst[2];
585 double min[2], max[2];
587 for (
unsigned int il = 0; il < 2; il++) {
590 std::string title = Form(
"final calibration const dist (%s): %s; entaRS (#alpha); entries",
m_label[il].data(),
m_runExp.data());
592 hnewconst[il]->SetTitle(Form(
"%s", title.data()));
594 title = Form(
"old calibration const dist (%s): %s; entaRS (#alpha); entries",
m_label[il].data(),
m_runExp.data());
596 holdconst[il]->SetTitle(Form(
"%s", title.data()));
598 for (
unsigned int iea = 0; iea < nbins; iea++) {
600 holdconst[il]->SetBinContent(iea + 1, prev);
601 hnewconst[il]->SetBinContent(iea + 1,
m_onedcors[il][iea]);
603 min[il] = hnewconst[il]->GetMinimum();
604 max[il] = hnewconst[il]->GetMaximum();
608 if (max[1] < max[0])max[1] = max[0];
609 if (min[1] > min[0])min[1] = min[0];
611 gStyle->SetOptStat(
"ne");
612 TCanvas cfconst(
"cfconst",
"Final calirbation constants", 800, 400);
613 cfconst.Divide(2, 1);
615 for (
int il = 0; il < 2; il++) {
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");
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");
629 cfconst.SaveAs(Form(
"cdcdedx_1dcell_fconsts%s.pdf",
m_suffix.data()));
630 cfconst.SaveAs(Form(
"cdcdedx_1dcell_fconsts%s.root",
m_suffix.data()));
632 for (
int il = 0; il < 2; il++) {
633 delete hnewconst[il];
634 delete holdconst[il];
642 TCanvas cstats(
"cstats",
"cstats", 1000, 500);
643 cstats.SetBatch(kTRUE);
647 auto hestats = getObjectPtr<TH1I>(
"hestats");
649 hestats->SetName(Form(
"hestats_%s",
m_suffix.data()));
650 hestats->SetStats(0);
651 hestats->DrawCopy(
"");
655 auto htstats = getObjectPtr<TH1I>(
"htstats");
657 hestats->SetName(Form(
"htstats_%s",
m_suffix.data()));
658 htstats->DrawCopy(
"");
659 hestats->SetStats(0);
661 cstats.Print(Form(
"cdcdedx_1dcell_stats_%s.pdf",
m_suffix.data()));
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
std::vector< int > m_eaBinLocal
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
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.