9#include <cdc/calibration/CDCdEdx/CDCDedxWireGainAlgorithm.h>
11#include <cdc/geometry/CDCGeometryPar.h>
37 setDescription(
"A calibration algorithm for CDC dE/dx wire gains");
50 B2FATAL(
"There is no valid payload for BadWires and/or Wirgain");
56 vector<int>* wire = 0;
57 ttree->SetBranchAddress(
"wire", &wire);
59 vector<double>* dedxhit = 0;
60 ttree->SetBranchAddress(
"dedxhit", &dedxhit);
63 ttree->SetBranchAddress(
"costh", &costh);
66 map<int, vector<double>> wirededx;
69 array<TH1D*, 2> hdedxL;
70 string label[2] = {
"IL",
"OL"};
72 for (
int il = 0; il < 2; il++)
75 int slWireBoundary = 2240;
77 for (
int i = 0; i < ttree->GetEntries(); ++i) {
79 if (costh > -0.55 && costh < 0.82) {
80 for (
unsigned int j = 0; j < wire->size(); ++j) {
81 int jwire = wire->at(j);
82 double jhitdedx = dedxhit->at(j);
83 wirededx[jwire].push_back(jhitdedx);
85 if (jwire < slWireBoundary) hdedxL[0]->Fill(jhitdedx);
86 else hdedxL[1]->Fill(jhitdedx);
93 for (
unsigned int jw = 0; jw < c_nwireCDC; ++jw)
94 if (wirededx[jw].size() <= 100) minstat++;
99 array<unsigned int, 2> minbinL, maxbinL;
100 for (
int il = 0; il < 2; il++) {
102 hdedxL[il]->SetTitle(Form(
"%s(%s);%d;%d", label[il].data(),
m_suffix.data(), minbinL[il], maxbinL[il]));
105 vector<double> vrel_mean;
106 vector<double> vdedx_means;
108 vector<TH1D*> hdedxhit(c_nwireCDC);
110 B2INFO(
"Creating CDCGeometryPar object");
113 vector<double> layermean(c_maxNSenseLayers);
115 int activelayers = 0;
116 double layeravg = 0.0;
119 for (
unsigned int il = 0; il < c_maxNSenseLayers; ++il) {
124 for (
unsigned int iw = 0; iw < cdcgeo.
nWiresInLayer(il); ++iw) {
129 for (
unsigned int ih = 0; ih < wirededx[jwire].size(); ++ih) {
130 hdedxhit[jwire]->Fill(wirededx[jwire][ih]);
133 unsigned int minbin, maxbin;
137 if (jwire < slWireBoundary) {
145 hdedxhit[jwire]->SetTitle(Form(
"dedxhit-dist, wire: %d (%s);%d;%d", jwire,
m_suffix.data(), minbin, maxbin));
148 if (
m_DBBadWires->getBadWireStatus(jwire) == kTRUE) dedxmean = 0.0;
151 vrel_mean.push_back(dedxmean);
155 vdedx_means.push_back(dedxmean * prewg);
156 B2INFO(
"merged-wireGain: [" << jwire <<
"], prewgvious = " << prewg <<
", rel = " << dedxmean <<
", merged = " << vdedx_means.at(
158 }
else vdedx_means.push_back(dedxmean);
161 if (vdedx_means.at(jwire) > 0) {
162 layermean[il] += vdedx_means.at(jwire);
167 if (activewires > 0) layermean[il] /= activewires;
168 else layermean[il] = 1.0;
171 unsigned int layerBoundary = 14;
172 if (il >= layerBoundary && layermean[il] > 0) {
173 layeravg += layermean[il];
180 if (activelayers > 0) layeravg /= activelayers;
182 for (
unsigned int iw = 0; iw < c_nwireCDC; ++iw) {
183 vrel_mean.at(iw) /= layeravg;
184 vdedx_means.at(iw) /= layeravg;
219 if (cruns == 0) B2INFO(
"CDCDedxWireGain: start exp " << expRun.first <<
" and run " << expRun.second <<
"");
224 int estart = erStart.first;
225 int rstart = erStart.second;
228 int rend = erEnd.second;
235 else m_suffix = Form(
"e%d_r%dr%d", estart, rstart, rend);
242 B2INFO(
"dE/dx Calibration done for " << vdedx_means.size() <<
" CDC wires");
252 if (hdedxhit->Integral() < 100)
return 1.0;
254 if (binlow <= 0 || binhigh > hdedxhit->GetNbinsX())
return 1.0;
256 double binweights = 0., sumofbc = 0.;
257 for (
int ibin = binlow; ibin <= binhigh; ibin++) {
258 double bcdedx = hdedxhit->GetBinContent(ibin);
260 binweights += (bcdedx * hdedxhit->GetBinCenter(ibin));
264 if (sumofbc > 0)
return binweights / sumofbc;
273 double sum = hdedxhit->Integral();
274 if (sum <= 0 || hdedxhit->GetNbinsX() <= 0) {
275 binlow = 1; binhigh = 1;
279 binlow = 1.0; binhigh = 1.0;
280 double sumPer5 = 0.0, sumPer75 = 0.0;
281 for (
int ibin = 1; ibin <= hdedxhit->GetNbinsX(); ibin++) {
282 double bcdedx = hdedxhit->GetBinContent(ibin);
299 TCanvas cldedx(
"cldedx",
"IL/OL dedxhit dist", 900, 400);
302 for (
int il = 0; il < 2; il++) {
305 int minbin = stoi(hdedxL[il]->GetXaxis()->GetTitle());
306 int maxbin = stoi(hdedxL[il]->GetYaxis()->GetTitle());
307 double lowedge = hdedxL[il]->GetXaxis()->GetBinLowEdge(minbin);
308 double upedge = hdedxL[il]->GetXaxis()->GetBinUpEdge(maxbin);
310 hdedxL[il]->SetFillColor(kYellow);
311 hdedxL[il]->SetTitle(Form(
"%s, trunc(%0.02f - %0.02f);dedxhit;entries", hdedxL[il]->GetTitle(), lowedge, upedge));
312 hdedxL[il]->Draw(
"histo");
314 TH1D* hdedxLC = (TH1D*)hdedxL[il]->Clone(Form(
"%s_c", hdedxL[il]->GetName()));
315 hdedxLC->GetXaxis()->SetRange(minbin, maxbin);
316 hdedxLC->SetFillColor(kAzure + 1);
317 hdedxLC->Draw(
"same histo");
320 cldedx.SaveAs(Form(
"cdcdedx_wgcal_layerdedx_%s.pdf",
m_suffix.data()));
327 TCanvas ctmp(Form(
"cdcdedx_%s",
m_suffix.data()),
"", 1200, 1200);
329 ctmp.SetBatch(kTRUE);
332 psname << Form(
"cdcdedx_wgcal_%s.pdf[",
m_suffix.data());
333 ctmp.Print(psname.str().c_str());
335 psname << Form(
"cdcdedx_wgcal_%s.pdf",
m_suffix.data());
337 for (
unsigned int iw = 0; iw < hist.size(); iw++) {
339 int minbin = stoi(hist[iw]->GetXaxis()->GetTitle());
340 int maxbin = stoi(hist[iw]->GetYaxis()->GetTitle());
342 hist[iw]->SetFillColor(kYellow - 9);
343 hist[iw]->SetTitle(Form(
"%s, rel. #mu_{trunc} %0.03f;dedxhit;entries", hist[iw]->GetTitle(), vrel_mean.at(iw)));
346 hist[iw]->SetLineColor(kRed);
347 hist[iw]->SetLineWidth(2);
349 ctmp.cd(iw % 16 + 1);
350 hist[iw]->DrawCopy();
352 TH1D* hdedxhitC = (TH1D*)hist[iw]->Clone(Form(
"%sC", hist[iw]->GetName()));
353 hdedxhitC->GetXaxis()->SetRange(minbin, maxbin);
354 hdedxhitC->SetFillColor(kAzure + 1);
355 hdedxhitC->DrawCopy(
"same histo");
357 if (((iw + 1) % 16 == 0) || iw == (hist.size() - 1)) {
358 ctmp.Print(psname.str().c_str());
366 psname << Form(
"cdcdedx_wgcal_%s.pdf]",
m_suffix.data());
367 ctmp.Print(psname.str().c_str());
375 TCanvas cwconst(
"cwconst",
"", 900, 500);
376 TCanvas cwconstvar(
"cwconstvar",
"", 500, 400);
378 array<TH1D*, 2> hconstpw, hconstpwvar;
380 for (
int i = 0; i < 2; i++) {
382 hconstpw[i] =
new TH1D(Form(
"hconstpw_%d_%s", i,
m_suffix.data()),
"", c_nwireCDC, -0.5, 14335.5);
383 hconstpw[i]->SetTitle(Form(
"wiregain const (%s); wire numbers;<dedxhit>",
m_suffix.data()));
385 && i == 0) hconstpw[i]->SetTitle(Form(
"merged wiregain rel-const (%s), avg = %0.03f; wire numbers;<dedxhit>",
m_suffix.data(),
388 hconstpwvar[i] =
new TH1D(Form(
"hconstpwvar_%s",
m_suffix.data()),
"", 400, -0.5, 2.5);
389 hconstpwvar[i]->SetTitle(Form(
"wiregain const (%s); wire gains; nentries",
m_suffix.data()));
391 && i == 0) hconstpwvar[i]->SetTitle(Form(
"merged wiregain rel-const (%s), avg = %0.03f; wire gains; nentries",
m_suffix.data(),
394 for (
unsigned int iw = 0; iw < c_nwireCDC; iw++) {
396 double gain = vdedx_means.at(iw);
397 if (
m_isMerge && i == 1) gain = vrel_mean.at(iw);
398 hconstpw[i]->SetBinContent(iw + 1, gain);
399 hconstpwvar[i]->Fill(gain);
401 if (iw % 500 == 0) hconstpw[i]->GetXaxis()->SetBinLabel(iw + 1, Form(
"w%d", iw + 1));
404 hconstpw[i]->SetLineColor(i * 2 + 2);
405 hconstpw[i]->LabelsOption(
"u",
"X");
406 hconstpw[i]->GetYaxis()->SetRangeUser(-0.1, 3.5);
407 hconstpw[i]->LabelsDeflate();
409 hconstpwvar[i]->SetFillColor(i * 2 + 2);
410 hconstpwvar[i]->Scale(1 / hconstpwvar[i]->GetMaximum());
415 hconstpw[0]->Draw(
"");
416 if (
m_isMerge) hconstpw[1]->Draw(
"same");
419 hconstpwvar[0]->Draw(
"hist");
420 if (
m_isMerge) hconstpwvar[1]->Draw(
"hist same");
422 cwconst.SaveAs(Form(
"cdcdedx_wgcal_wireconst_%s.pdf",
m_suffix.data()));
423 cwconstvar.SaveAs(Form(
"cdcdedx_wgcal_wireconstvar_%s.pdf",
m_suffix.data()));
430 TH1D hlayeravg(Form(
"hlayeravg_%s",
m_suffix.data()),
"", layermean.size(), -0.5, 55.5);
431 hlayeravg.SetTitle(Form(
"layer gain avg (%s); layer numbers;<dedxhit>",
m_suffix.data()));
433 for (
unsigned int il = 0; il < layermean.size(); il++) {
434 hlayeravg.SetBinContent(il + 1, layermean[il]);
435 if (il % 2 == 0 || il == layermean.size() - 1) hlayeravg.GetXaxis()->SetBinLabel(il + 1, Form(
"L%d", il));
438 TCanvas clayeravg(
"clayeravg",
"clayeravg", 800, 500);
439 clayeravg.SetGridy(1);
441 gStyle->SetOptStat(
"ne");
442 hlayeravg.LabelsOption(
"u",
"X");
443 hlayeravg.SetLineColor(kBlue);
444 hlayeravg.GetYaxis()->SetRangeUser(-0.1, 3.5);
445 hlayeravg.SetTitle(Form(
"%s, avg = %0.03f (abs)", hlayeravg.GetTitle(), layeravg));
446 hlayeravg.LabelsDeflate();
448 TLine* tl =
new TLine(-0.5, layeravg, 55.5, layeravg);
449 tl->SetLineColor(kRed);
450 tl->DrawClone(
"same");
451 clayeravg.SaveAs(Form(
"cdcdedx_wgcal_layeravg_%s.pdf",
m_suffix.data()));
461 TCanvas clconst(
"clconst",
"", 800, 500);
462 clconst.Divide(2, 2);
463 clconst.SetBatch(kTRUE);
465 stringstream psnameL;
466 psnameL << Form(
"cdcdedx_wgcal_layerconst_%s.pdf[",
m_suffix.data());
467 clconst.Print(psnameL.str().c_str());
468 psnameL.str(
""); psnameL << Form(
"cdcdedx_wgcal_layerconst_%s.pdf",
m_suffix.data());
472 for (
unsigned int il = 0; il < layermean.size(); ++il) {
475 TH1D hconstpl(Form(
"hconstpwvar_l%d_%s", il,
m_suffix.data()),
"", nwires, jwire, jwire + nwires);
476 hconstpl.SetTitle(Form(
"abs-const, layer: %d (%s); wire numbers;<dedxhit>", il,
m_suffix.data()));
478 for (
unsigned int iw = 0; iw < nwires; ++iw) {
479 hconstpl.SetBinContent(iw + 1, vdedx_means.at(jwire));
481 if (iw % 10 == 0) hconstpl.GetXaxis()->SetBinLabel(iw + 1, Form(
"w%d", jwire));
483 if (iw % 15 == 0) hconstpl.GetXaxis()->SetBinLabel(iw + 1, Form(
"w%d", jwire));
488 double lmean = layermean.at(il) / layeravg;
490 clconst.cd(il % 4 + 1);
491 gStyle->SetOptStat(
"ne");
493 hconstpl.SetTitle(Form(
"%s, avg = %0.03f", hconstpl.GetTitle(), lmean));
495 if (il < 8) hconstpl.GetYaxis()->SetRangeUser(-0.1, 4.0);
496 else hconstpl.GetYaxis()->SetRangeUser(-0.1, 2.0);
497 hconstpl.SetFillColor(kAzure - 1);
498 hconstpl.LabelsOption(
"u",
"X");
499 hconstpl.DrawCopy(
"hist");
501 TLine* tlc =
new TLine(jwire - nwires, lmean, jwire, lmean);
502 tlc->SetLineColor(kRed);
503 tlc->DrawClone(
"same");
505 if ((il + 1) % 4 == 0 || (il + 1) == layermean.size()) {
506 clconst.Print(psnameL.str().c_str());
514 psnameL << Form(
"cdcdedx_wgcal_layerconst_%s.pdf]",
m_suffix.data());
515 clconst.Print(psnameL.str().c_str());
522 TCanvas cstats(
"cstats",
"cstats", 800, 400);
523 cstats.SetBatch(kTRUE);
529 hestats->SetName(Form(
"hestats_%s",
m_suffix.data()));
530 hestats->SetStats(0);
531 hestats->DrawCopy(
"");
537 htstats->SetName(Form(
"htstats_%s",
m_suffix.data()));
538 htstats->SetStats(0);
539 htstats->DrawCopy(
"");
542 cstats.Print(Form(
"cdcdedx_wgcal_stats_%s.pdf",
m_suffix.data()));
void plotWGPerLayer(const std::vector< double > &vdedx_means, const std::vector< double > &layermean, double layeravg)
function to draw WG per layer
void plotLayerGain(const std::vector< double > &layermean, double layeravg)
function to draw layer gains
void plotLayerDist(std::array< TH1D *, 2 > hdedxL)
function to draw dE/dx for inner/outer layer
bool m_isMerge
merge payload at the time of calibration
double m_truncMax
max trunc range for mean
double m_truncMin
min trunc range for mean
void getTruncatedBins(TH1D *hdedxhit, unsigned int &binlow, unsigned int &binhigh)
function to get bins of truncation from histogram
void plotWireGain(const std::vector< double > &vdedx_means, const std::vector< double > &vrel_mean, double layeravg)
function to draw wire gains
void getExpRunInfo()
function to get extract calibration run/exp
void createPayload(const std::vector< double > &vdedx_tmeans)
function to finally store new payload after full calibration
void plotWireDist(const std::vector< TH1D * > &hist, const std::vector< double > &vrel_mean)
function to draw dE/dx histograms for each wire
bool m_isMakePlots
Save arithmetic and truncated mean for the 'dedx' values.
DBObjPtr< CDCDedxBadWires > m_DBBadWires
Bad wire DB object.
std::string m_suffix
suffix string to separate plots
double getTruncationMean(TH1D *hdedxhit, int binlow, int binhigh)
function to get mean of truncation from histogram
int m_dedxBins
number of bins for dedx histogram
DBObjPtr< CDCGeometry > m_cdcGeo
Geometry of CDC.
CDCDedxWireGainAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
int m_exp
exp no to set SL boundaries
DBObjPtr< CDCDedxWireGain > m_DBWireGains
Wire gain DB object.
void plotEventStats()
function to draw statstics
virtual EResult calibrate() override
Wire gain algorithm.
double m_dedxMax
max dedx range for wiregain cal
double m_dedxMin
min dedx range for wiregain cal
bool m_isWireTruc
method of trunc range for mean
dE/dx wire gain calibration constants
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.
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.