9 #include <reconstruction/calibration/CDCDedxBadWireAlgorithm.h>
20 c_nwireCDC(c_nSenseWires),
46 B2FATAL(
"There is no valid payload for BadWire and/or Wirgain");
49 auto ttree = getObjectPtr<TTree>(
"tree");
52 vector<int>* wire = 0;
53 ttree->SetBranchAddress(
"wire", &wire);
55 vector<double>* hitvar = 0;
56 if (
m_isADC) ttree->SetBranchAddress(
"adccorr", &hitvar);
57 else ttree->SetBranchAddress(
"dedxhit", &hitvar);
63 hvarall.SetTitle(Form(
"dist %s; %s; %s",
m_suffix.data(),
m_varName.data(),
"entries"));
65 map<int, vector<double>> vhitvar;
67 for (
int i = 0; i < ttree->GetEntries(); ++i) {
69 for (
unsigned int ih = 0; ih < wire->size(); ++ih) {
70 double ivalue = hitvar->at(ih);
71 vhitvar[wire->at(ih)].push_back(ivalue);
81 for (
unsigned int jw = 0; jw <
c_nwireCDC; ++jw)
82 if (vhitvar[jw].size() <= 100) minstat++;
86 map<int, vector<double>> qapars;
87 vector<double> vdefectwires, vbadwires, vdeadwires;
89 for (
unsigned int jw = 0; jw <
c_nwireCDC; ++jw) {
90 int ncount = 0, tcount = 0;
92 for (
unsigned int jh = 0; jh < vhitvar[jw].size(); ++jh) {
93 double jvalue = vhitvar[jw][jh];
100 bool badwire =
false;
102 qapars[0].push_back(0);
103 qapars[1].push_back(0);
104 qapars[2].push_back(0);
107 nmean = nmean / ncount;
111 for (
unsigned int kh = 0; kh < vhitvar[jw].size(); ++kh) {
112 double kvalue = vhitvar[jw][kh];
113 if (kvalue <
m_varMax) nrms += pow(kvalue - nmean, 2);
116 nrms =
sqrt(nrms / ncount);
119 double badfrac = 0.0;
120 if (tcount > 0) badfrac = (1.0 * tcount) / (tcount + ncount);
123 qapars[0].push_back(nmean);
124 qapars[1].push_back(nrms);
125 qapars[2].push_back(badfrac);
129 vdefectwires.push_back(0.0);
130 if (ncount == 0) vdeadwires.push_back(jw);
131 else vbadwires.push_back(jw);
132 }
else vdefectwires.push_back(1.0);
151 B2INFO(
"dE/dx Badwire Calibration done: " << vdefectwires.size() <<
" wires");
166 if (cruns == 0) B2INFO(
"CDCDedxBadWires: start exp " << expRun.first <<
" and run " << expRun.second <<
"");
171 int estart = erStart.first;
172 int rstart = erStart.second;
175 int rend = erEnd.second;
180 else m_suffix = Form(
"e%d_r%dr%d", estart, rstart, rend);
186 map<
int, vector<double>>& vhitvar)
189 TList* bdlist =
new TList();
190 bdlist->SetName(
"badwires");
192 TList* goodlist =
new TList();
193 goodlist->SetName(
"goodwires");
195 TList* hflist =
new TList();
196 hflist->SetName(
"highfracwires");
198 for (
unsigned int jw = 0; jw <
c_nwireCDC; ++jw) {
203 hvarhf->SetTitle(Form(
"%s, wire = %d; %s; entries",
m_suffix.data(), jw,
m_varName.data()));
205 int ncount = 0, tcount = 0;
207 for (
unsigned int jh = 0; jh < vhitvar[jw].size(); ++jh) {
208 double jvalue = vhitvar[jw][jh];
214 if (jvalue <
m_varMax * 10.) hvarhf->Fill(jvalue / 10.);
218 double badfrac = 0.0;
219 if (tcount > 0) badfrac = (1.0 * tcount) / (tcount + ncount);
220 hvar->SetTitle(Form(
"%s, wire = %d; %s; %0.01f",
m_suffix.data(), jw,
m_varName.data(), badfrac * 100));
223 if (count(inwires.begin(), inwires.end(), jw)) isbad =
true;
227 hvar->SetLineWidth(2);
228 hvar->SetLineColor(kRed);
235 if (hvar->Integral() > 100) goodlist->Add(hvar);
251 string listname = list->GetName();
252 string sfx = Form(
"%s_%s", listname.data(),
m_suffix.data());
254 TCanvas ctmp(Form(
"cdcdedx_%s", sfx.data()),
"", 1200, 1200);
256 ctmp.SetBatch(kTRUE);
259 psname << Form(
"cdcdedx_bdcal_%s.pdf[", sfx.data());
260 ctmp.Print(psname.str().c_str());
262 psname << Form(
"cdcdedx_bdcal_%s.pdf", sfx.data());
264 for (
int ih = 0; ih < list->GetSize(); ih++) {
266 TH1D* hist = (TH1D*)list->At(ih);
268 double frac = stod(hist->GetYaxis()->GetTitle());
270 TPaveText* pinfo =
new TPaveText(0.40, 0.63, 0.89, 0.89,
"NBNDC");
273 pinfo->AddText(Form(
"N: %0.00f", hist->Integral()));
274 pinfo->AddText(Form(
"hf: %0.00f%%(%0.00f%%)", frac,
m_fracThres * 100));
277 ctmp.cd(ih % 16 + 1);
278 hist->GetYaxis()->SetTitle(
"entries");
279 hist->SetFillColor(color);
284 if (listname ==
"badwires") {
285 TH1D* histhf = (TH1D*)hflist->At(ih);
286 if (hist->GetMaximum() < histhf->GetMaximum()) hist->SetMaximum(histhf->GetMaximum() * 1.05);
287 histhf->SetFillColor(kGray);
289 histhf->Draw(
"same");
292 if (((ih + 1) % 16 == 0) || ih == (list->GetSize() - 1)) {
293 ctmp.Print(psname.str().c_str());
299 psname << Form(
"cdcdedx_bdcal_%s.pdf]", sfx.data());
300 ctmp.Print(psname.str().c_str());
307 TCanvas cmap(Form(
"cmap_%s",
m_suffix.data()),
"", 800, 800);
308 cmap.SetTitle(
"CDC dE/dx bad wire status");
312 hxyAll->SetTitle(Form(
"wire status map (%s)",
m_suffix.data()));
320 hxyBad->Draw(
"same");
327 hxyDead->Draw(
"same");
330 int ndefect = nbad + ndead;
331 auto leg =
new TLegend(0.68, 0.80, 0.90, 0.92);
332 leg->SetBorderSize(0);
333 leg->SetLineWidth(3);
334 leg->SetHeader(Form(
"total defective: %d (~%0.02f%%)", ndefect, 100.*(ndefect) /
c_nwireCDC));
335 leg->AddEntry(hxyBad, Form(
"bad #rightarrow %d", nbad),
"p");
336 leg->AddEntry(hxyDead, Form(
"dead #rightarrow %d", ndead),
"p");
339 gStyle->SetLegendTextSize(0.025);
340 TPaveText* pt =
new TPaveText(-0.30, -1.47, -0.31, -1.30,
"br");
343 TText* text = pt->AddText(
"CDC-wire map: counter-clockwise and start from +x");
344 text->SetTextColor(kGray + 1);
347 cmap.SaveAs(Form(
"cdcdedx_bdcal_wiremap_%s.pdf",
m_suffix.data()));
358 B2INFO(
"Creating CDCGeometryPar object");
361 TH2F* temp =
new TH2F(Form(
"temp_%s_%s",
m_suffix.data(), suffix.data()),
"", 2400, -1.2, 1.2, 2400, -1.2, 1.2);
365 for (
unsigned int ilay = 0; ilay < c_maxNSenseLayers; ++ilay) {
366 for (
unsigned int iwire = 0; iwire < cdcgeo.
nWiresInLayer(ilay); ++iwire) {
368 double phi = 2.*TMath::Pi() * (float(iwire) / float(cdcgeo.
nWiresInLayer(ilay)));
369 double radius = cdcgeo.
senseWireR(ilay) / 100.;
370 double x = radius * cos(phi);
371 double y = radius * sin(phi);
372 if (suffix ==
"all") {
376 if (count(inwires.begin(), inwires.end(), jwire)) {
390 string qaname[3] = {
"mean",
"rms",
"high_fraction"};
395 for (
int iqa = 0; iqa < 3; iqa++) {
397 TH1D histqa(Form(
"%s_%s", qaname[iqa].data(),
m_suffix.data()),
"",
c_nwireCDC, -0.5, 14335.5);
399 for (
unsigned int jw = 0; jw <
c_nwireCDC; jw++) {
400 if (iqa == 2) histqa.SetBinContent(jw + 1, qapars[iqa][jw] * 100);
401 else histqa.SetBinContent(jw + 1, qapars[iqa][jw]);
404 TCanvas c_pars(Form(
"c_pars_%d", iqa),
"", 800, 600);
408 histqa.SetTitle(Form(
"%s vs wires (%s); wire ; %s", qaname[iqa].data(),
m_suffix.data(), qaname[iqa].data()));
412 TLine* lmin =
new TLine(-0.5, linemin[iqa], 14335.5, linemin[iqa]);
413 lmin->SetLineColor(kRed);
415 TLine* lmax =
new TLine(-0.5, linemax[iqa], 14335.5, linemax[iqa]);
416 lmax->SetLineColor(kRed);
419 c_pars.Print(Form(
"cdcdedx_bdcal_%s_%s.root", qaname[iqa].data(),
m_suffix.data()));
420 c_pars.Print(Form(
"cdcdedx_bdcal_%s_%s.pdf", qaname[iqa].data(),
m_suffix.data()));
431 TCanvas cstats(
"cstats",
"cstats", 1000, 500);
432 cstats.SetBatch(kTRUE);
436 auto hestats = getObjectPtr<TH1I>(
"hestats");
438 hestats->SetName(Form(
"htstats_%s",
m_suffix.data()));
439 hestats->SetStats(0);
440 hestats->DrawCopy(
"");
444 auto htstats = getObjectPtr<TH1I>(
"htstats");
446 htstats->SetName(Form(
"htstats_%s",
m_suffix.data()));
447 htstats->SetStats(0);
448 htstats->DrawCopy(
"");
451 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
void plotWireDist(const std::vector< double > &inwires, std::map< int, std::vector< double >> &vhitvar)
function to draw per wire plots
double m_varMin
min range for input variable
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
double m_arms
average rms of dedx for all wires
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.
void plotQaPars(std::map< int, std::vector< double >> &qapars)
function to plot the QA (decision) parameters
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 plotEventStats()
function to draw the stats
virtual EResult calibrate() override
cdcdedx badwire algorithm
double m_amean
average mean of dedx for all 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.
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)
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
double sqrt(double a)
sqrt for double
Abstract base class for different kinds of events.