Belle II Software development
CDCDedxBadWireAlgorithm.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <cdc/calibration/CDCdEdx/CDCDedxBadWireAlgorithm.h>
10
11#include <cdc/geometry/CDCGeometryPar.h>
12
13#include <TCanvas.h>
14#include <TLegend.h>
15#include <TLine.h>
16#include <TMath.h>
17#include <TStyle.h>
18
19using namespace Belle2;
20using namespace CDC;
21using namespace std;
22
23//-----------------------------------------------------------------
24// Implementation
25//-----------------------------------------------------------------
27 CalibrationAlgorithm("CDCDedxElectronCollector"),
28 c_nwireCDC(c_nSenseWires),
29 m_isMakePlots(true),
30 m_isADC(false),
31 m_varBins(100),
32 m_varMin(0.0),
33 m_varMax(7.0),
34 m_meanThres(1.0),
35 m_rmsThres(1.0),
36 m_fracThres(1.0),
37 m_varName("hitdedx"),
38 m_suffix("")
39{
40 // Set module properties
41 setDescription("A calibration algorithm for CDC dE/dx bad wires");
42}
43
44//-----------------------------------------------------------------
45// Run the calibration
46//-----------------------------------------------------------------
48{
49
51
52 //old wg for book-keeping previous bad+dead
53 if (!m_DBBadWires.isValid() || !m_DBWireGains.isValid())
54 B2FATAL("There is no valid payload for BadWire and/or Wirgain");
55
56 // Get data objects
57 auto ttree = getObjectPtr<TTree>("tree");
58 if (ttree->GetEntries() < 1000) return c_NotEnoughData;
59
60 vector<int>* wire = 0;
61 ttree->SetBranchAddress("wire", &wire);
62
63 vector<double>* hitvar = 0;
64 if (m_isADC) ttree->SetBranchAddress("adccorr", &hitvar);
65 else ttree->SetBranchAddress("dedxhit", &hitvar);
66
67 if (m_isADC) m_varName = "hitadc";
68 m_suffix = Form("%s_%s", m_varName.data(), m_suffix.data());
69
70 //IL = Inner Layer and OL = Outer Layer
71 array<TH1D*, 2> hvarL;
72 string label[2] = {"IL", "OL"};
73 for (int il = 0; il < 2; il++) {
74 hvarL[il] = new TH1D(Form("hvar%s_%s", label[il].data(), m_suffix.data()), "", m_varBins, m_varMin, m_varMax);
75 hvarL[il]->SetTitle(Form("dist (%s) %s; %s; %s", label[il].data(), m_suffix.data(), m_varName.data(), "entries"));
76 }
77
78 map<int, vector<double>> vhitvar;
79
80 m_slWireBoundary = (m_exp < 40) ? 1280 : 2240;
81
82 for (int i = 0; i < ttree->GetEntries(); ++i) {
83 ttree->GetEvent(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);
88
89 if (jwire < m_slWireBoundary) hvarL[0]->Fill(ivalue);
90 else hvarL[1]->Fill(ivalue);
91 }
92 }
93
94 m_amean_IL = hvarL[0]->GetMean();
95 m_arms_IL = hvarL[0]->GetRMS();
96 m_amean_OL = hvarL[1]->GetMean();
97 m_arms_OL = hvarL[1]->GetRMS();
98
99 // Commenting minstat cut on March 2024,
100 //it is not worth to skip the calibration if there is more than 5% bad/dead wires
101
102 //return if >5% bad wire or null histogram
103 // int minstat = 0;
104 // for (unsigned int jw = 0; jw < c_nwireCDC; ++jw)
105 // if (vhitvar[jw].size() <= 100) minstat++;
106 // if (minstat > 0.05 * c_nwireCDC) return c_NotEnoughData;
107
108 if (m_amean_IL == 0 || m_arms_IL == 0 || m_amean_OL == 0 || m_arms_OL == 0) return c_NotEnoughData;
109
110 map<int, vector<double>> qapars;
111 vector<double> vdefectwires, vbadwires, vdeadwires;
112
113 for (unsigned int jw = 0; jw < c_nwireCDC; ++jw) {
114 int ncount = 0, tcount = 0;
115 double nmean = 0.;
116 for (unsigned int jh = 0; jh < vhitvar[jw].size(); ++jh) {
117 double jvalue = vhitvar[jw][jh];
118 if (jvalue < m_varMax) {
119 ncount++;
120 nmean += jvalue;
121 } else tcount++;
122 }
123
124 bool badwire = false;
125 if (ncount < 100) {
126 qapars[0].push_back(0);
127 qapars[1].push_back(0);
128 qapars[2].push_back(0);
129 badwire = true; //partial dead
130 } else {
131 nmean = nmean / ncount;
132 if (int(jw) < m_slWireBoundary) {
133 if (abs(nmean - m_amean_IL) / m_amean_IL > m_meanThres) badwire = true;
134 } else
135 {if (abs(nmean - m_amean_OL) / m_amean_OL > m_meanThres) badwire = true;}
136
137 double nrms = 0.;
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);
141 }
142
143 nrms = sqrt(nrms / ncount);
144 if (int(jw) < m_slWireBoundary)
145 {if (abs(nrms - m_arms_IL) / m_arms_IL > m_rmsThres) badwire = true;}
146 else
147 {if (abs(nrms - m_arms_OL) / m_arms_OL > m_rmsThres) badwire = true;}
148
149 double badfrac = 0.0;
150 if (tcount > 0) badfrac = (1.0 * tcount) / (tcount + ncount);
151 if (badfrac > m_fracThres) badwire = true;
152
153 qapars[0].push_back(nmean);
154 qapars[1].push_back(nrms);
155 qapars[2].push_back(badfrac);
156 }
157
158 if (badwire) {
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);
163 }
164
165
166 if (m_isMakePlots) {
167 //1. plot bad and good wire plots.
168 plotWireDist(vbadwires, vhitvar);
169
170 //2. plots wire status map
171 plotBadWireMap(vbadwires, vdeadwires);
172
173 //3. plot control parameters histograms
174 plotQaPars(qapars);
175
176 //4. plot statistics related histograms
178 }
179
180 // Save payloads
181 B2INFO("dE/dx Badwire Calibration done: " << vdefectwires.size() << " wires");
182 CDCDedxBadWires* c_badwires = new CDCDedxBadWires(vdefectwires);
183 saveCalibration(c_badwires, "CDCDedxBadWires");
184
185 m_suffix.clear();
186
187 return c_OK;
188}
189
190//------------------------------------
192{
193
194 int cruns = 0;
195 for (auto expRun : getRunList()) {
196 if (cruns == 0) B2INFO("CDCDedxBadWires: start exp " << expRun.first << " and run " << expRun.second << "");
197 cruns++;
198 }
199
200 const auto erStart = getRunList()[0];
201 int estart = erStart.first;
202 int rstart = erStart.second;
203
204 const auto erEnd = getRunList()[cruns - 1];
205 int rend = erEnd.second;
206
207 m_exp = estart;
208
209 updateDBObjPtrs(1, rstart, estart);
210
211 if (m_suffix.length() > 0) m_suffix = Form("%s_e%d_r%dr%d", m_suffix.data(), estart, rstart, rend);
212 else m_suffix = Form("e%d_r%dr%d", estart, rstart, rend);
213}
214
215
216//------------------------------------
217void CDCDedxBadWireAlgorithm::plotWireDist(const vector<double>& inwires,
218 map<int, vector<double>>& vhitvar)
219{
220
221 TList* bdlist = new TList();
222 bdlist->SetName("badwires");
223
224 TList* goodlist = new TList();
225 goodlist->SetName("goodwires");
226
227 TList* hflist = new TList();
228 hflist->SetName("highfracwires");
229
230 for (unsigned int jw = 0; jw < c_nwireCDC; ++jw) {
231
232 TH1D* hvar = new TH1D(Form("%s_wire%d", m_suffix.data(), jw), "", m_varBins, m_varMin, m_varMax);
233 hvar->SetUniqueID(jw);
234
235 TH1D* hvarhf = new TH1D(Form("hf%s_wire%d", m_suffix.data(), jw), "", m_varBins, m_varMin, m_varMax);
236 hvarhf->SetUniqueID(jw);
237 hvarhf->SetTitle(Form("%s, wire = %d; %s; entries", m_suffix.data(), jw, m_varName.data()));
238
239 int ncount = 0, tcount = 0;
240
241 for (unsigned int jh = 0; jh < vhitvar[jw].size(); ++jh) {
242 double jvalue = vhitvar[jw][jh];
243 if (jvalue < m_varMax) {
244 ncount++;
245 hvar->Fill(jvalue);
246 } else {
247 tcount++;
248 if (jvalue < m_varMax * 10.) hvarhf->Fill(jvalue / 10.);
249 }
250 }
251
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));
255
256 bool isbad = false;
257 if (count(inwires.begin(), inwires.end(), jw)) isbad = true;
258
259 double oldwg = m_DBWireGains->getWireGain(jw);
260 if (oldwg == 0) {
261 hvar->SetLineWidth(2);
262 hvar->SetLineColor(kRed);
263 }
264
265 if (isbad) {
266 bdlist->Add(hvar);
267 hflist->Add(hvarhf);
268 } else {
269 if (hvar->Integral() > 100) goodlist->Add(hvar);
270 }
271 }
272
273 printCanvas(bdlist, hflist, kYellow - 9);
274 printCanvas(goodlist, hflist, kGreen);
275
276 delete bdlist;
277 delete goodlist;
278 delete hflist;
279}
280
281//------------------------------------
282void CDCDedxBadWireAlgorithm::printCanvas(TList* list, TList* hflist, Color_t color)
283{
284
285 string listname = list->GetName();
286 string sfx = Form("%s_%s", listname.data(), m_suffix.data());
287
288 TCanvas ctmp(Form("cdcdedx_%s", sfx.data()), "", 1200, 1200);
289 ctmp.Divide(4, 4);
290 ctmp.SetBatch(kTRUE);
291
292 stringstream psname;
293 psname << Form("cdcdedx_bdcal_%s.pdf[", sfx.data());
294 ctmp.Print(psname.str().c_str());
295 psname.str("");
296 psname << Form("cdcdedx_bdcal_%s.pdf", sfx.data());
297
298 for (int ih = 0; ih < list->GetSize(); ih++) {
299
300 TH1D* hist = (TH1D*)list->At(ih);
301 int jw = hist->GetUniqueID();
302
303 double frac = stod(hist->GetYaxis()->GetTitle());
304
305 double amean = (jw < m_slWireBoundary) ? m_amean_IL : m_amean_OL;
306 double arms = (jw < m_slWireBoundary) ? m_arms_IL : m_arms_OL;
307
308 TPaveText* pinfo = new TPaveText(0.40, 0.63, 0.89, 0.89, "NBNDC");
309
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));
314 setTextCosmetics(pinfo, 0.04258064);
315
316 ctmp.cd(ih % 16 + 1);
317 hist->GetYaxis()->SetTitle("entries");
318 hist->SetFillColor(color);
319 hist->SetStats(0);
320 hist->Draw();
321 pinfo->Draw("same");
322
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);
327 histhf->SetStats(0);
328 histhf->Draw("same");
329 }
330
331 if (((ih + 1) % 16 == 0) || ih == (list->GetSize() - 1)) {
332 ctmp.Print(psname.str().c_str());
333 ctmp.Clear("D");
334 }
335 }
336
337 psname.str("");
338 psname << Form("cdcdedx_bdcal_%s.pdf]", sfx.data());
339 ctmp.Print(psname.str().c_str());
340}
341
342//------------------------------------
343void CDCDedxBadWireAlgorithm::plotBadWireMap(const vector<double>& vbadwires, const vector<double>& vdeadwires)
344{
345
346 TCanvas cmap(Form("cmap_%s", m_suffix.data()), "", 800, 800);
347 cmap.SetTitle("CDC dE/dx bad wire status");
348
349 int total = 0;
350 TH2F* hxyAll = getHistoPattern(vbadwires, "all", total);
351 hxyAll->SetTitle(Form("wire status map (%s)", m_suffix.data()));
352 setHistCosmetics(hxyAll, kGray);
353 hxyAll->Draw();
354
355 int nbad = 0.0;
356 TH2F* hxyBad = getHistoPattern(vbadwires, "bad", nbad);
357 if (hxyBad) {
358 setHistCosmetics(hxyBad, kRed);
359 hxyBad->Draw("same");
360 }
361
362 int ndead = 0.0;
363 TH2F* hxyDead = getHistoPattern(vdeadwires, "dead", ndead);
364 if (hxyDead) {
365 setHistCosmetics(hxyDead, kBlack);
366 hxyDead->Draw("same");
367 }
368
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");
376 leg->Draw();
377
378 gStyle->SetLegendTextSize(0.025);
379 TPaveText* pt = new TPaveText(-0.30, -1.47, -0.31, -1.30, "br");
380 setTextCosmetics(pt, 0.02258064);
381
382 TText* text = pt->AddText("CDC-wire map: counter-clockwise and start from +x");
383 text->SetTextColor(kGray + 1);
384 pt->Draw("same");
385
386 cmap.SaveAs(Form("cdcdedx_bdcal_wiremap_%s.pdf", m_suffix.data()));
387
388 delete hxyAll;
389 delete hxyBad;
390 delete hxyDead;
391}
392
393//------------------------------------
394TH2F* CDCDedxBadWireAlgorithm::getHistoPattern(const vector<double>& inwires, const string& suffix, int& total)
395{
396
397 B2INFO("Creating CDCGeometryPar object");
399
400 TH2F* temp = new TH2F(Form("temp_%s_%s", m_suffix.data(), suffix.data()), "", 2400, -1.2, 1.2, 2400, -1.2, 1.2);
401
402 int jwire = -1;
403 total = 0;
404 for (unsigned int ilay = 0; ilay < c_maxNSenseLayers; ++ilay) {
405 for (unsigned int iwire = 0; iwire < cdcgeo.nWiresInLayer(ilay); ++iwire) {
406 jwire++;
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") {
412 total++;
413 temp->Fill(x, y);
414 } else {
415 if (count(inwires.begin(), inwires.end(), jwire)) {
416 temp->Fill(x, y);
417 total++;
418 }
419 }
420 }
421 }
422 return temp;
423}
424
425//------------------------------------
426void CDCDedxBadWireAlgorithm::plotQaPars(map<int, vector<double>>& qapars)
427{
428
429 string qaname[3] = {"mean", "rms", "high_fraction"};
430
431 double lineminIL[3] = {m_amean_IL* (1 - m_meanThres), m_arms_IL* (1 - m_rmsThres), m_fracThres * 100};
432 double linemaxIL[3] = {m_amean_IL* (1 + m_meanThres), m_arms_IL* (1 + m_rmsThres), m_fracThres * 100};
433
434 double lineminOL[3] = {m_amean_OL* (1 - m_meanThres), m_arms_OL* (1 - m_rmsThres), m_fracThres * 100};
435 double linemaxOL[3] = {m_amean_OL* (1 + m_meanThres), m_arms_OL* (1 + m_rmsThres), m_fracThres * 100};
436
437 for (int iqa = 0; iqa < 3; iqa++) {
438
439 TH1D histqa(Form("%s_%s", qaname[iqa].data(), m_suffix.data()), "", c_nwireCDC, -0.5, 14335.5);
440
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]);
444 }
445
446 TCanvas c_pars(Form("c_pars_%d", iqa), "", 800, 600);
447 c_pars.cd();
448 gPad->SetGridy();
449
450 histqa.SetTitle(Form("%s vs wires (%s); wire ; %s", qaname[iqa].data(), m_suffix.data(), qaname[iqa].data()));
451 histqa.SetStats(0);
452 histqa.Draw();
453
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");
466
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()));
469
470 delete lmaxIL;
471 delete lminIL;
472 }
473}
474
475//------------------------------------
477{
478
479 TCanvas cstats("cstats", "cstats", 1000, 500);
480 cstats.SetBatch(kTRUE);
481 cstats.Divide(2, 1);
482
483 cstats.cd(1);
484 auto hestats = getObjectPtr<TH1I>("hestats");
485 if (hestats) {
486 hestats->SetName(Form("htstats_%s", m_suffix.data()));
487 hestats->SetStats(0);
488 hestats->DrawCopy("");
489 }
490
491 cstats.cd(2);
492 auto htstats = getObjectPtr<TH1I>("htstats");
493 if (htstats) {
494 htstats->SetName(Form("htstats_%s", m_suffix.data()));
495 htstats->SetStats(0);
496 htstats->DrawCopy("");
497 }
498
499 cstats.Print(Form("cdcdedx_bdcal_qastats_%s.pdf", m_suffix.data()));
500}
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
Definition beamHelpers.h:28
Abstract base class for different kinds of events.
STL namespace.