9#include <cdc/calibration/CDCdEdx/CDCDedxBadWireAlgorithm.h>
11#include <cdc/geometry/CDCGeometryPar.h>
54 B2FATAL(
"There is no valid payload for BadWire and/or Wirgain");
60 vector<int>* wire = 0;
61 ttree->SetBranchAddress(
"wire", &wire);
63 vector<double>* hitvar = 0;
64 if (
m_isADC) ttree->SetBranchAddress(
"adccorr", &hitvar);
65 else ttree->SetBranchAddress(
"dedxhit", &hitvar);
71 array<TH1D*, 2> hvarL;
72 string label[2] = {
"IL",
"OL"};
73 for (
int il = 0; il < 2; il++) {
75 hvarL[il]->SetTitle(Form(
"dist (%s) %s; %s; %s", label[il].data(),
m_suffix.data(),
m_varName.data(),
"entries"));
78 map<int, vector<double>> vhitvar;
82 for (
int i = 0; i < ttree->GetEntries(); ++i) {
84 for (
unsigned int ih = 0; ih < wire->size(); ++ih) {
85 int jwire = wire->at(ih);
86 double ivalue = hitvar->at(ih);
87 vhitvar[wire->at(ih)].push_back(ivalue);
90 else hvarL[1]->Fill(ivalue);
110 map<int, vector<double>> qapars;
111 vector<double> vdefectwires, vbadwires, vdeadwires;
113 for (
unsigned int jw = 0; jw <
c_nwireCDC; ++jw) {
114 int ncount = 0, tcount = 0;
116 for (
unsigned int jh = 0; jh < vhitvar[jw].size(); ++jh) {
117 double jvalue = vhitvar[jw][jh];
124 bool badwire =
false;
126 qapars[0].push_back(0);
127 qapars[1].push_back(0);
128 qapars[2].push_back(0);
131 nmean = nmean / ncount;
138 for (
unsigned int kh = 0; kh < vhitvar[jw].size(); ++kh) {
139 double kvalue = vhitvar[jw][kh];
140 if (kvalue <
m_varMax) nrms += pow(kvalue - nmean, 2);
143 nrms =
sqrt(nrms / ncount);
149 double badfrac = 0.0;
150 if (tcount > 0) badfrac = (1.0 * tcount) / (tcount + ncount);
153 qapars[0].push_back(nmean);
154 qapars[1].push_back(nrms);
155 qapars[2].push_back(badfrac);
159 vdefectwires.push_back(0.0);
160 if (ncount == 0) vdeadwires.push_back(jw);
161 else vbadwires.push_back(jw);
162 }
else vdefectwires.push_back(1.0);
181 B2INFO(
"dE/dx Badwire Calibration done: " << vdefectwires.size() <<
" wires");
196 if (cruns == 0) B2INFO(
"CDCDedxBadWires: start exp " << expRun.first <<
" and run " << expRun.second <<
"");
201 int estart = erStart.first;
202 int rstart = erStart.second;
205 int rend = erEnd.second;
212 else m_suffix = Form(
"e%d_r%dr%d", estart, rstart, rend);
218 map<
int, vector<double>>& vhitvar)
221 TList* bdlist =
new TList();
222 bdlist->SetName(
"badwires");
224 TList* goodlist =
new TList();
225 goodlist->SetName(
"goodwires");
227 TList* hflist =
new TList();
228 hflist->SetName(
"highfracwires");
230 for (
unsigned int jw = 0; jw <
c_nwireCDC; ++jw) {
233 hvar->SetUniqueID(jw);
236 hvarhf->SetUniqueID(jw);
237 hvarhf->SetTitle(Form(
"%s, wire = %d; %s; entries",
m_suffix.data(), jw,
m_varName.data()));
239 int ncount = 0, tcount = 0;
241 for (
unsigned int jh = 0; jh < vhitvar[jw].size(); ++jh) {
242 double jvalue = vhitvar[jw][jh];
248 if (jvalue <
m_varMax * 10.) hvarhf->Fill(jvalue / 10.);
252 double badfrac = 0.0;
253 if (tcount > 0) badfrac = (1.0 * tcount) / (tcount + ncount);
254 hvar->SetTitle(Form(
"%s, wire = %d; %s; %0.01f",
m_suffix.data(), jw,
m_varName.data(), badfrac * 100));
257 if (count(inwires.begin(), inwires.end(), jw)) isbad =
true;
261 hvar->SetLineWidth(2);
262 hvar->SetLineColor(kRed);
269 if (hvar->Integral() > 100) goodlist->Add(hvar);
285 string listname = list->GetName();
286 string sfx = Form(
"%s_%s", listname.data(),
m_suffix.data());
288 TCanvas ctmp(Form(
"cdcdedx_%s", sfx.data()),
"", 1200, 1200);
290 ctmp.SetBatch(kTRUE);
293 psname << Form(
"cdcdedx_bdcal_%s.pdf[", sfx.data());
294 ctmp.Print(psname.str().c_str());
296 psname << Form(
"cdcdedx_bdcal_%s.pdf", sfx.data());
298 for (
int ih = 0; ih < list->GetSize(); ih++) {
300 TH1D* hist = (TH1D*)list->At(ih);
301 int jw = hist->GetUniqueID();
303 double frac = stod(hist->GetYaxis()->GetTitle());
308 TPaveText* pinfo =
new TPaveText(0.40, 0.63, 0.89, 0.89,
"NBNDC");
310 pinfo->AddText(Form(
"#mu: %0.2f(%0.2f#pm%0.2f)", hist->GetMean(), amean,
m_meanThres * amean));
311 pinfo->AddText(Form(
"#sigma: %0.2f(%0.2f#pm%0.2f)", hist->GetRMS(), arms,
m_rmsThres * arms));
312 pinfo->AddText(Form(
"N: %0.00f", hist->Integral()));
313 pinfo->AddText(Form(
"hf: %0.00f%%(%0.00f%%)", frac,
m_fracThres * 100));
316 ctmp.cd(ih % 16 + 1);
317 hist->GetYaxis()->SetTitle(
"entries");
318 hist->SetFillColor(color);
323 if (listname ==
"badwires") {
324 TH1D* histhf = (TH1D*)hflist->At(ih);
325 if (hist->GetMaximum() < histhf->GetMaximum()) hist->SetMaximum(histhf->GetMaximum() * 1.05);
326 histhf->SetFillColor(kGray);
328 histhf->Draw(
"same");
331 if (((ih + 1) % 16 == 0) || ih == (list->GetSize() - 1)) {
332 ctmp.Print(psname.str().c_str());
338 psname << Form(
"cdcdedx_bdcal_%s.pdf]", sfx.data());
339 ctmp.Print(psname.str().c_str());
346 TCanvas cmap(Form(
"cmap_%s",
m_suffix.data()),
"", 800, 800);
347 cmap.SetTitle(
"CDC dE/dx bad wire status");
351 hxyAll->SetTitle(Form(
"wire status map (%s)",
m_suffix.data()));
359 hxyBad->Draw(
"same");
366 hxyDead->Draw(
"same");
369 int ndefect = nbad + ndead;
370 auto leg =
new TLegend(0.68, 0.80, 0.90, 0.92);
371 leg->SetBorderSize(0);
372 leg->SetLineWidth(3);
373 leg->SetHeader(Form(
"total defective: %d (~%0.02f%%)", ndefect, 100.*(ndefect) /
c_nwireCDC));
374 leg->AddEntry(hxyBad, Form(
"bad #rightarrow %d", nbad),
"p");
375 leg->AddEntry(hxyDead, Form(
"dead #rightarrow %d", ndead),
"p");
378 gStyle->SetLegendTextSize(0.025);
379 TPaveText* pt =
new TPaveText(-0.30, -1.47, -0.31, -1.30,
"br");
382 TText* text = pt->AddText(
"CDC-wire map: counter-clockwise and start from +x");
383 text->SetTextColor(kGray + 1);
386 cmap.SaveAs(Form(
"cdcdedx_bdcal_wiremap_%s.pdf",
m_suffix.data()));
397 B2INFO(
"Creating CDCGeometryPar object");
400 TH2F* temp =
new TH2F(Form(
"temp_%s_%s",
m_suffix.data(), suffix.data()),
"", 2400, -1.2, 1.2, 2400, -1.2, 1.2);
404 for (
unsigned int ilay = 0; ilay < c_maxNSenseLayers; ++ilay) {
405 for (
unsigned int iwire = 0; iwire < cdcgeo.
nWiresInLayer(ilay); ++iwire) {
407 double phi = 2.*TMath::Pi() * (float(iwire) / float(cdcgeo.
nWiresInLayer(ilay)));
408 double radius = cdcgeo.
senseWireR(ilay) / 100.;
409 double x = radius * cos(phi);
410 double y = radius * sin(phi);
411 if (suffix ==
"all") {
415 if (count(inwires.begin(), inwires.end(), jwire)) {
429 string qaname[3] = {
"mean",
"rms",
"high_fraction"};
437 for (
int iqa = 0; iqa < 3; iqa++) {
439 TH1D histqa(Form(
"%s_%s", qaname[iqa].data(),
m_suffix.data()),
"",
c_nwireCDC, -0.5, 14335.5);
441 for (
unsigned int jw = 0; jw <
c_nwireCDC; jw++) {
442 if (iqa == 2) histqa.SetBinContent(jw + 1, qapars[iqa][jw] * 100);
443 else histqa.SetBinContent(jw + 1, qapars[iqa][jw]);
446 TCanvas c_pars(Form(
"c_pars_%d", iqa),
"", 800, 600);
450 histqa.SetTitle(Form(
"%s vs wires (%s); wire ; %s", qaname[iqa].data(),
m_suffix.data(), qaname[iqa].data()));
454 TLine* lminIL =
new TLine(-0.5, lineminIL[iqa], 2239.5, lineminIL[iqa]);
455 lminIL->SetLineColor(kRed);
456 lminIL->Draw(
"same");
457 TLine* lmaxIL =
new TLine(-0.5, linemaxIL[iqa], 2239.5, linemaxIL[iqa]);
458 lmaxIL->SetLineColor(kRed);
459 lmaxIL->Draw(
"same");
460 TLine* lminOL =
new TLine(2239.5, lineminOL[iqa], 14335.5, lineminOL[iqa]);
461 lminOL->SetLineColor(kRed);
462 lminOL->Draw(
"same");
463 TLine* lmaxOL =
new TLine(2239.5, linemaxOL[iqa], 14335.5, linemaxOL[iqa]);
464 lmaxOL->SetLineColor(kRed);
465 lmaxOL->Draw(
"same");
467 c_pars.Print(Form(
"cdcdedx_bdcal_%s_%s.root", qaname[iqa].data(),
m_suffix.data()));
468 c_pars.Print(Form(
"cdcdedx_bdcal_%s_%s.pdf", qaname[iqa].data(),
m_suffix.data()));
479 TCanvas cstats(
"cstats",
"cstats", 1000, 500);
480 cstats.SetBatch(kTRUE);
486 hestats->SetName(Form(
"htstats_%s",
m_suffix.data()));
487 hestats->SetStats(0);
488 hestats->DrawCopy(
"");
494 htstats->SetName(Form(
"htstats_%s",
m_suffix.data()));
495 htstats->SetStats(0);
496 htstats->DrawCopy(
"");
499 cstats.Print(Form(
"cdcdedx_bdcal_qastats_%s.pdf",
m_suffix.data()));
void setHistCosmetics(TH2F *hist, Color_t color)
function to change histogram styles
void plotBadWireMap(const std::vector< double > &vbadwires, const std::vector< double > &vdeadwires)
function to plot wire status map (all, bad and dead)
double m_varMax
max range for input variable
double m_rmsThres
rms Threshold accepted for good wire
double m_varMin
min range for input variable
int m_slWireBoundary
Boundary between inner layers: SL0 (<40), SL0+SL1 (>=40)
double m_amean_IL
average mean of dedx for inner wires
std::string m_varName
std::string to set var name (adc or dedx)
void getExpRunInfo()
function to get extract calibration run/exp
unsigned int c_nwireCDC
number of wires in CDC
bool m_isMakePlots
produce plots for status
DBObjPtr< CDCDedxBadWires > m_DBBadWires
Badwire DB object.
std::string m_suffix
suffix std::string for naming plots
void setTextCosmetics(TPaveText *pt, double size)
function to change text styles
CDCDedxBadWireAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
DBObjPtr< CDCGeometry > m_cdcGeo
Geometry of CDC.
int m_exp
exp no to set SL boundaries
double m_meanThres
mean Threshold accepted for good wire
DBObjPtr< CDCDedxWireGain > m_DBWireGains
Wiregain DB object.
double m_fracThres
high-frac Threshold accepted for good wire
int m_varBins
number of bins for input variable
void printCanvas(TList *list, TList *hflist, Color_t color)
function to print canvas
void plotQaPars(std::map< int, std::vector< double > > &qapars)
function to plot the QA (decision) parameters
void plotEventStats()
function to draw the stats
virtual EResult calibrate() override
cdcdedx badwire algorithm
void plotWireDist(const std::vector< double > &inwires, std::map< int, std::vector< double > > &vhitvar)
function to draw per wire plots
double m_arms_IL
average rms of dedx for inner wires
double m_arms_OL
average rms of dedx for outer wires
double m_amean_OL
average mean of dedx for outer wires
bool m_isADC
Use adc if(true) else dedx for calibration.
TH2F * getHistoPattern(const std::vector< double > &inwires, const std::string &suffix, int &total)
function to get wire map with input file (all, bad and dead)
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.
double senseWireR(int layerId) const
Returns radius of sense wire in each layer.
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...
double sqrt(double a)
sqrt for double
Abstract base class for different kinds of events.