Belle II Software  release-08-01-10
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 <reconstruction/calibration/CDCDedxBadWireAlgorithm.h>
10 
11 using namespace Belle2;
12 using namespace CDC;
13 using namespace std;
14 
15 //-----------------------------------------------------------------
16 // Implementation
17 //-----------------------------------------------------------------
19  CalibrationAlgorithm("CDCDedxElectronCollector"),
20  c_nwireCDC(c_nSenseWires),
21  m_isMakePlots(true),
22  m_isADC(false),
23  m_varBins(100),
24  m_varMin(0.0),
25  m_varMax(7.0),
26  m_meanThres(1.0),
27  m_rmsThres(1.0),
28  m_fracThres(1.0),
29  m_varName("hitdedx"),
30  m_suffix("")
31 {
32  // Set module properties
33  setDescription("A calibration algorithm for CDC dE/dx bad wires");
34 }
35 
36 //-----------------------------------------------------------------
37 // Run the calibration
38 //-----------------------------------------------------------------
40 {
41 
42  getExpRunInfo();
43 
44  //old wg for book-keeping previous bad+dead
45  if (!m_DBBadWires.isValid() || !m_DBWireGains.isValid())
46  B2FATAL("There is no valid payload for BadWire and/or Wirgain");
47 
48  // Get data objects
49  auto ttree = getObjectPtr<TTree>("tree");
50  if (ttree->GetEntries() < 1000) return c_NotEnoughData;
51 
52  vector<int>* wire = 0;
53  ttree->SetBranchAddress("wire", &wire);
54 
55  vector<double>* hitvar = 0;
56  if (m_isADC) ttree->SetBranchAddress("adccorr", &hitvar);
57  else ttree->SetBranchAddress("dedxhit", &hitvar);
58 
59  if (m_isADC) m_varName = "hitadc";
60  m_suffix = Form("%s_%s", m_varName.data(), m_suffix.data());
61 
62  TH1D hvarall(Form("hvarall_%s", m_suffix.data()), "", m_varBins, m_varMin, m_varMax);
63  hvarall.SetTitle(Form("dist %s; %s; %s", m_suffix.data(), m_varName.data(), "entries"));
64 
65  map<int, vector<double>> vhitvar;
66 
67  for (int i = 0; i < ttree->GetEntries(); ++i) {
68  ttree->GetEvent(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);
72  hvarall.Fill(ivalue);
73  }
74  }
75 
76  m_amean = hvarall.GetMean();
77  m_arms = hvarall.GetRMS();
78 
79  // Commenting minstat cut on March 2024,
80  //it is not worth to skip the calibration if there is more than 5% bad/dead wires
81 
82  //return if >5% bad wire or null histogram
83  // int minstat = 0;
84  // for (unsigned int jw = 0; jw < c_nwireCDC; ++jw)
85  // if (vhitvar[jw].size() <= 100) minstat++;
86  // if (minstat > 0.05 * c_nwireCDC) return c_NotEnoughData;
87 
88  if (m_amean == 0 || m_arms == 0) return c_NotEnoughData;
89 
90  map<int, vector<double>> qapars;
91  vector<double> vdefectwires, vbadwires, vdeadwires;
92 
93  for (unsigned int jw = 0; jw < c_nwireCDC; ++jw) {
94  int ncount = 0, tcount = 0;
95  double nmean = 0.;
96  for (unsigned int jh = 0; jh < vhitvar[jw].size(); ++jh) {
97  double jvalue = vhitvar[jw][jh];
98  if (jvalue < m_varMax) {
99  ncount++;
100  nmean += jvalue;
101  } else tcount++;
102  }
103 
104  bool badwire = false;
105  if (ncount < 100) {
106  qapars[0].push_back(0);
107  qapars[1].push_back(0);
108  qapars[2].push_back(0);
109  badwire = true; //partial dead
110  } else {
111  nmean = nmean / ncount;
112  if (abs(nmean - m_amean) / m_amean > m_meanThres) badwire = true;
113 
114  double nrms = 0.;
115  for (unsigned int kh = 0; kh < vhitvar[jw].size(); ++kh) {
116  double kvalue = vhitvar[jw][kh];
117  if (kvalue < m_varMax) nrms += pow(kvalue - nmean, 2);
118  }
119 
120  nrms = sqrt(nrms / ncount);
121  if (abs(nrms - m_arms) / m_arms > m_rmsThres) badwire = true;
122 
123  double badfrac = 0.0;
124  if (tcount > 0) badfrac = (1.0 * tcount) / (tcount + ncount);
125  if (badfrac > m_fracThres) badwire = true;
126 
127  qapars[0].push_back(nmean);
128  qapars[1].push_back(nrms);
129  qapars[2].push_back(badfrac);
130  }
131 
132  if (badwire) {
133  vdefectwires.push_back(0.0);
134  if (ncount == 0) vdeadwires.push_back(jw);
135  else vbadwires.push_back(jw);
136  } else vdefectwires.push_back(1.0);
137  }
138 
139 
140  if (m_isMakePlots) {
141  //1. plot bad and good wire plots.
142  plotWireDist(vbadwires, vhitvar);
143 
144  //2. plots wire status map
145  plotBadWireMap(vbadwires, vdeadwires);
146 
147  //3. plot control parameters histograms
148  plotQaPars(qapars);
149 
150  //4. plot statistics related histograms
151  plotEventStats();
152  }
153 
154  // Save payloads
155  B2INFO("dE/dx Badwire Calibration done: " << vdefectwires.size() << " wires");
156  CDCDedxBadWires* c_badwires = new CDCDedxBadWires(vdefectwires);
157  saveCalibration(c_badwires, "CDCDedxBadWires");
158 
159  m_suffix.clear();
160 
161  return c_OK;
162 }
163 
164 //------------------------------------
166 {
167 
168  int cruns = 0;
169  for (auto expRun : getRunList()) {
170  if (cruns == 0) B2INFO("CDCDedxBadWires: start exp " << expRun.first << " and run " << expRun.second << "");
171  cruns++;
172  }
173 
174  const auto erStart = getRunList()[0];
175  int estart = erStart.first;
176  int rstart = erStart.second;
177 
178  const auto erEnd = getRunList()[cruns - 1];
179  int rend = erEnd.second;
180 
181  updateDBObjPtrs(1, rstart, estart);
182 
183  if (m_suffix.length() > 0) m_suffix = Form("%s_e%d_r%dr%d", m_suffix.data(), estart, rstart, rend);
184  else m_suffix = Form("e%d_r%dr%d", estart, rstart, rend);
185 }
186 
187 
188 //------------------------------------
189 void CDCDedxBadWireAlgorithm::plotWireDist(const vector<double>& inwires,
190  map<int, vector<double>>& vhitvar)
191 {
192 
193  TList* bdlist = new TList();
194  bdlist->SetName("badwires");
195 
196  TList* goodlist = new TList();
197  goodlist->SetName("goodwires");
198 
199  TList* hflist = new TList();
200  hflist->SetName("highfracwires");
201 
202  for (unsigned int jw = 0; jw < c_nwireCDC; ++jw) {
203 
204  TH1D* hvar = new TH1D(Form("%s_wire%d", m_suffix.data(), jw), "", m_varBins, m_varMin, m_varMax);
205 
206  TH1D* hvarhf = new TH1D(Form("hf%s_wire%d", m_suffix.data(), jw), "", m_varBins, m_varMin, m_varMax);
207  hvarhf->SetTitle(Form("%s, wire = %d; %s; entries", m_suffix.data(), jw, m_varName.data()));
208 
209  int ncount = 0, tcount = 0;
210 
211  for (unsigned int jh = 0; jh < vhitvar[jw].size(); ++jh) {
212  double jvalue = vhitvar[jw][jh];
213  if (jvalue < m_varMax) {
214  ncount++;
215  hvar->Fill(jvalue);
216  } else {
217  tcount++;
218  if (jvalue < m_varMax * 10.) hvarhf->Fill(jvalue / 10.);
219  }
220  }
221 
222  double badfrac = 0.0;
223  if (tcount > 0) badfrac = (1.0 * tcount) / (tcount + ncount);
224  hvar->SetTitle(Form("%s, wire = %d; %s; %0.01f", m_suffix.data(), jw, m_varName.data(), badfrac * 100));
225 
226  bool isbad = false;
227  if (count(inwires.begin(), inwires.end(), jw)) isbad = true;
228 
229  double oldwg = m_DBWireGains->getWireGain(jw);
230  if (oldwg == 0) {
231  hvar->SetLineWidth(2);
232  hvar->SetLineColor(kRed);
233  }
234 
235  if (isbad) {
236  bdlist->Add(hvar);
237  hflist->Add(hvarhf);
238  } else {
239  if (hvar->Integral() > 100) goodlist->Add(hvar);
240  }
241  }
242 
243  printCanvas(bdlist, hflist, kYellow - 9);
244  printCanvas(goodlist, hflist, kGreen);
245 
246  delete bdlist;
247  delete goodlist;
248  delete hflist;
249 }
250 
251 //------------------------------------
252 void CDCDedxBadWireAlgorithm::printCanvas(TList* list, TList* hflist, Color_t color)
253 {
254 
255  string listname = list->GetName();
256  string sfx = Form("%s_%s", listname.data(), m_suffix.data());
257 
258  TCanvas ctmp(Form("cdcdedx_%s", sfx.data()), "", 1200, 1200);
259  ctmp.Divide(4, 4);
260  ctmp.SetBatch(kTRUE);
261 
262  stringstream psname;
263  psname << Form("cdcdedx_bdcal_%s.pdf[", sfx.data());
264  ctmp.Print(psname.str().c_str());
265  psname.str("");
266  psname << Form("cdcdedx_bdcal_%s.pdf", sfx.data());
267 
268  for (int ih = 0; ih < list->GetSize(); ih++) {
269 
270  TH1D* hist = (TH1D*)list->At(ih);
271 
272  double frac = stod(hist->GetYaxis()->GetTitle());
273 
274  TPaveText* pinfo = new TPaveText(0.40, 0.63, 0.89, 0.89, "NBNDC");
275  pinfo->AddText(Form("#mu: %0.2f(%0.2f#pm%0.2f)", hist->GetMean(), m_amean, m_meanThres * m_amean));
276  pinfo->AddText(Form("#sigma: %0.2f(%0.2f#pm%0.2f)", hist->GetRMS(), m_arms, m_rmsThres * m_arms));
277  pinfo->AddText(Form("N: %0.00f", hist->Integral()));
278  pinfo->AddText(Form("hf: %0.00f%%(%0.00f%%)", frac, m_fracThres * 100));
279  setTextCosmetics(pinfo, 0.04258064);
280 
281  ctmp.cd(ih % 16 + 1);
282  hist->GetYaxis()->SetTitle("entries");
283  hist->SetFillColor(color);
284  hist->SetStats(0);
285  hist->Draw();
286  pinfo->Draw("same");
287 
288  if (listname == "badwires") {
289  TH1D* histhf = (TH1D*)hflist->At(ih);
290  if (hist->GetMaximum() < histhf->GetMaximum()) hist->SetMaximum(histhf->GetMaximum() * 1.05);
291  histhf->SetFillColor(kGray);
292  histhf->SetStats(0);
293  histhf->Draw("same");
294  }
295 
296  if (((ih + 1) % 16 == 0) || ih == (list->GetSize() - 1)) {
297  ctmp.Print(psname.str().c_str());
298  ctmp.Clear("D");
299  }
300  }
301 
302  psname.str("");
303  psname << Form("cdcdedx_bdcal_%s.pdf]", sfx.data());
304  ctmp.Print(psname.str().c_str());
305 }
306 
307 //------------------------------------
308 void CDCDedxBadWireAlgorithm::plotBadWireMap(const vector<double>& vbadwires, const vector<double>& vdeadwires)
309 {
310 
311  TCanvas cmap(Form("cmap_%s", m_suffix.data()), "", 800, 800);
312  cmap.SetTitle("CDC dE/dx bad wire status");
313 
314  int total = 0;
315  TH2F* hxyAll = getHistoPattern(vbadwires, "all", total);
316  hxyAll->SetTitle(Form("wire status map (%s)", m_suffix.data()));
317  setHistCosmetics(hxyAll, kGray);
318  hxyAll->Draw();
319 
320  int nbad = 0.0;
321  TH2F* hxyBad = getHistoPattern(vbadwires, "bad", nbad);
322  if (hxyBad) {
323  setHistCosmetics(hxyBad, kRed);
324  hxyBad->Draw("same");
325  }
326 
327  int ndead = 0.0;
328  TH2F* hxyDead = getHistoPattern(vdeadwires, "dead", ndead);
329  if (hxyDead) {
330  setHistCosmetics(hxyDead, kBlack);
331  hxyDead->Draw("same");
332  }
333 
334  int ndefect = nbad + ndead;
335  auto leg = new TLegend(0.68, 0.80, 0.90, 0.92);
336  leg->SetBorderSize(0);
337  leg->SetLineWidth(3);
338  leg->SetHeader(Form("total defective: %d (~%0.02f%%)", ndefect, 100.*(ndefect) / c_nwireCDC));
339  leg->AddEntry(hxyBad, Form("bad #rightarrow %d", nbad), "p");
340  leg->AddEntry(hxyDead, Form("dead #rightarrow %d", ndead), "p");
341  leg->Draw();
342 
343  gStyle->SetLegendTextSize(0.025);
344  TPaveText* pt = new TPaveText(-0.30, -1.47, -0.31, -1.30, "br");
345  setTextCosmetics(pt, 0.02258064);
346 
347  TText* text = pt->AddText("CDC-wire map: counter-clockwise and start from +x");
348  text->SetTextColor(kGray + 1);
349  pt->Draw("same");
350 
351  cmap.SaveAs(Form("cdcdedx_bdcal_wiremap_%s.pdf", m_suffix.data()));
352 
353  delete hxyAll;
354  delete hxyBad;
355  delete hxyDead;
356 }
357 
358 //------------------------------------
359 TH2F* CDCDedxBadWireAlgorithm::getHistoPattern(const vector<double>& inwires, const string& suffix, int& total)
360 {
361 
362  B2INFO("Creating CDCGeometryPar object");
364 
365  TH2F* temp = new TH2F(Form("temp_%s_%s", m_suffix.data(), suffix.data()), "", 2400, -1.2, 1.2, 2400, -1.2, 1.2);
366 
367  int jwire = -1;
368  total = 0;
369  for (unsigned int ilay = 0; ilay < c_maxNSenseLayers; ++ilay) {
370  for (unsigned int iwire = 0; iwire < cdcgeo.nWiresInLayer(ilay); ++iwire) {
371  jwire++;
372  double phi = 2.*TMath::Pi() * (float(iwire) / float(cdcgeo.nWiresInLayer(ilay)));
373  double radius = cdcgeo.senseWireR(ilay) / 100.;
374  double x = radius * cos(phi);
375  double y = radius * sin(phi);
376  if (suffix == "all") {
377  total++;
378  temp->Fill(x, y);
379  } else {
380  if (count(inwires.begin(), inwires.end(), jwire)) {
381  temp->Fill(x, y);
382  total++;
383  }
384  }
385  }
386  }
387  return temp;
388 }
389 
390 //------------------------------------
391 void CDCDedxBadWireAlgorithm::plotQaPars(map<int, vector<double>>& qapars)
392 {
393 
394  string qaname[3] = {"mean", "rms", "high_fraction"};
395 
396  double linemin[3] = {m_amean* (1 - m_meanThres), m_arms* (1 - m_rmsThres), m_fracThres * 100};
397  double linemax[3] = {m_amean* (1 + m_meanThres), m_arms* (1 + m_rmsThres), m_fracThres * 100};
398 
399  for (int iqa = 0; iqa < 3; iqa++) {
400 
401  TH1D histqa(Form("%s_%s", qaname[iqa].data(), m_suffix.data()), "", c_nwireCDC, -0.5, 14335.5);
402 
403  for (unsigned int jw = 0; jw < c_nwireCDC; jw++) {
404  if (iqa == 2) histqa.SetBinContent(jw + 1, qapars[iqa][jw] * 100);
405  else histqa.SetBinContent(jw + 1, qapars[iqa][jw]);
406  }
407 
408  TCanvas c_pars(Form("c_pars_%d", iqa), "", 800, 600);
409  c_pars.cd();
410  gPad->SetGridy();
411 
412  histqa.SetTitle(Form("%s vs wires (%s); wire ; %s", qaname[iqa].data(), m_suffix.data(), qaname[iqa].data()));
413  histqa.SetStats(0);
414  histqa.Draw();
415 
416  TLine* lmin = new TLine(-0.5, linemin[iqa], 14335.5, linemin[iqa]);
417  lmin->SetLineColor(kRed);
418  lmin->Draw("same");
419  TLine* lmax = new TLine(-0.5, linemax[iqa], 14335.5, linemax[iqa]);
420  lmax->SetLineColor(kRed);
421  lmax->Draw("same");
422 
423  c_pars.Print(Form("cdcdedx_bdcal_%s_%s.root", qaname[iqa].data(), m_suffix.data()));
424  c_pars.Print(Form("cdcdedx_bdcal_%s_%s.pdf", qaname[iqa].data(), m_suffix.data()));
425 
426  delete lmax;
427  delete lmin;
428  }
429 }
430 
431 //------------------------------------
433 {
434 
435  TCanvas cstats("cstats", "cstats", 1000, 500);
436  cstats.SetBatch(kTRUE);
437  cstats.Divide(2, 1);
438 
439  cstats.cd(1);
440  auto hestats = getObjectPtr<TH1I>("hestats");
441  if (hestats) {
442  hestats->SetName(Form("htstats_%s", m_suffix.data()));
443  hestats->SetStats(0);
444  hestats->DrawCopy("");
445  }
446 
447  cstats.cd(2);
448  auto htstats = getObjectPtr<TH1I>("htstats");
449  if (htstats) {
450  htstats->SetName(Form("htstats_%s", m_suffix.data()));
451  htstats->SetStats(0);
452  htstats->DrawCopy("");
453  }
454 
455  cstats.Print(Form("cdcdedx_bdcal_qastats_%s.pdf", m_suffix.data()));
456 }
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
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.