Belle II Software  release-08-00-05
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  //return if >5% bad wire or null histogram
80  int minstat = 0;
81  for (unsigned int jw = 0; jw < c_nwireCDC; ++jw)
82  if (vhitvar[jw].size() <= 100) minstat++;
83 
84  if (minstat > 0.05 * c_nwireCDC || m_amean == 0 || m_arms == 0) return c_NotEnoughData;
85 
86  map<int, vector<double>> qapars;
87  vector<double> vdefectwires, vbadwires, vdeadwires;
88 
89  for (unsigned int jw = 0; jw < c_nwireCDC; ++jw) {
90  int ncount = 0, tcount = 0;
91  double nmean = 0.;
92  for (unsigned int jh = 0; jh < vhitvar[jw].size(); ++jh) {
93  double jvalue = vhitvar[jw][jh];
94  if (jvalue < m_varMax) {
95  ncount++;
96  nmean += jvalue;
97  } else tcount++;
98  }
99 
100  bool badwire = false;
101  if (ncount < 100) {
102  qapars[0].push_back(0);
103  qapars[1].push_back(0);
104  qapars[2].push_back(0);
105  badwire = true; //partial dead
106  } else {
107  nmean = nmean / ncount;
108  if (abs(nmean - m_amean) / m_amean > m_meanThres) badwire = true;
109 
110  double nrms = 0.;
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);
114  }
115 
116  nrms = sqrt(nrms / ncount);
117  if (abs(nrms - m_arms) / m_arms > m_rmsThres) badwire = true;
118 
119  double badfrac = 0.0;
120  if (tcount > 0) badfrac = (1.0 * tcount) / (tcount + ncount);
121  if (badfrac > m_fracThres) badwire = true;
122 
123  qapars[0].push_back(nmean);
124  qapars[1].push_back(nrms);
125  qapars[2].push_back(badfrac);
126  }
127 
128  if (badwire) {
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);
133  }
134 
135 
136  if (m_isMakePlots) {
137  //1. plot bad and good wire plots.
138  plotWireDist(vbadwires, vhitvar);
139 
140  //2. plots wire status map
141  plotBadWireMap(vbadwires, vdeadwires);
142 
143  //3. plot control parameters histograms
144  plotQaPars(qapars);
145 
146  //4. plot statistics related histograms
147  plotEventStats();
148  }
149 
150  // Save payloads
151  B2INFO("dE/dx Badwire Calibration done: " << vdefectwires.size() << " wires");
152  CDCDedxBadWires* c_badwires = new CDCDedxBadWires(vdefectwires);
153  saveCalibration(c_badwires, "CDCDedxBadWires");
154 
155  m_suffix.clear();
156 
157  return c_OK;
158 }
159 
160 //------------------------------------
162 {
163 
164  int cruns = 0;
165  for (auto expRun : getRunList()) {
166  if (cruns == 0) B2INFO("CDCDedxBadWires: start exp " << expRun.first << " and run " << expRun.second << "");
167  cruns++;
168  }
169 
170  const auto erStart = getRunList()[0];
171  int estart = erStart.first;
172  int rstart = erStart.second;
173 
174  const auto erEnd = getRunList()[cruns - 1];
175  int rend = erEnd.second;
176 
177  updateDBObjPtrs(1, rstart, estart);
178 
179  if (m_suffix.length() > 0) m_suffix = Form("%s_e%d_r%dr%d", m_suffix.data(), estart, rstart, rend);
180  else m_suffix = Form("e%d_r%dr%d", estart, rstart, rend);
181 }
182 
183 
184 //------------------------------------
185 void CDCDedxBadWireAlgorithm::plotWireDist(const vector<double>& inwires,
186  map<int, vector<double>>& vhitvar)
187 {
188 
189  TList* bdlist = new TList();
190  bdlist->SetName("badwires");
191 
192  TList* goodlist = new TList();
193  goodlist->SetName("goodwires");
194 
195  TList* hflist = new TList();
196  hflist->SetName("highfracwires");
197 
198  for (unsigned int jw = 0; jw < c_nwireCDC; ++jw) {
199 
200  TH1D* hvar = new TH1D(Form("%s_wire%d", m_suffix.data(), jw), "", m_varBins, m_varMin, m_varMax);
201 
202  TH1D* hvarhf = new TH1D(Form("hf%s_wire%d", m_suffix.data(), jw), "", m_varBins, m_varMin, m_varMax);
203  hvarhf->SetTitle(Form("%s, wire = %d; %s; entries", m_suffix.data(), jw, m_varName.data()));
204 
205  int ncount = 0, tcount = 0;
206 
207  for (unsigned int jh = 0; jh < vhitvar[jw].size(); ++jh) {
208  double jvalue = vhitvar[jw][jh];
209  if (jvalue < m_varMax) {
210  ncount++;
211  hvar->Fill(jvalue);
212  } else {
213  tcount++;
214  if (jvalue < m_varMax * 10.) hvarhf->Fill(jvalue / 10.);
215  }
216  }
217 
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));
221 
222  bool isbad = false;
223  if (count(inwires.begin(), inwires.end(), jw)) isbad = true;
224 
225  double oldwg = m_DBWireGains->getWireGain(jw);
226  if (oldwg == 0) {
227  hvar->SetLineWidth(2);
228  hvar->SetLineColor(kRed);
229  }
230 
231  if (isbad) {
232  bdlist->Add(hvar);
233  hflist->Add(hvarhf);
234  } else {
235  if (hvar->Integral() > 100) goodlist->Add(hvar);
236  }
237  }
238 
239  printCanvas(bdlist, hflist, kYellow - 9);
240  printCanvas(goodlist, hflist, kGreen);
241 
242  delete bdlist;
243  delete goodlist;
244  delete hflist;
245 }
246 
247 //------------------------------------
248 void CDCDedxBadWireAlgorithm::printCanvas(TList* list, TList* hflist, Color_t color)
249 {
250 
251  string listname = list->GetName();
252  string sfx = Form("%s_%s", listname.data(), m_suffix.data());
253 
254  TCanvas ctmp(Form("cdcdedx_%s", sfx.data()), "", 1200, 1200);
255  ctmp.Divide(4, 4);
256  ctmp.SetBatch(kTRUE);
257 
258  stringstream psname;
259  psname << Form("cdcdedx_bdcal_%s.pdf[", sfx.data());
260  ctmp.Print(psname.str().c_str());
261  psname.str("");
262  psname << Form("cdcdedx_bdcal_%s.pdf", sfx.data());
263 
264  for (int ih = 0; ih < list->GetSize(); ih++) {
265 
266  TH1D* hist = (TH1D*)list->At(ih);
267 
268  double frac = stod(hist->GetYaxis()->GetTitle());
269 
270  TPaveText* pinfo = new TPaveText(0.40, 0.63, 0.89, 0.89, "NBNDC");
271  pinfo->AddText(Form("#mu: %0.2f(%0.2f#pm%0.2f)", hist->GetMean(), m_amean, m_meanThres * m_amean));
272  pinfo->AddText(Form("#sigma: %0.2f(%0.2f#pm%0.2f)", hist->GetRMS(), m_arms, m_rmsThres * m_arms));
273  pinfo->AddText(Form("N: %0.00f", hist->Integral()));
274  pinfo->AddText(Form("hf: %0.00f%%(%0.00f%%)", frac, m_fracThres * 100));
275  setTextCosmetics(pinfo, 0.04258064);
276 
277  ctmp.cd(ih % 16 + 1);
278  hist->GetYaxis()->SetTitle("entries");
279  hist->SetFillColor(color);
280  hist->SetStats(0);
281  hist->Draw();
282  pinfo->Draw("same");
283 
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);
288  histhf->SetStats(0);
289  histhf->Draw("same");
290  }
291 
292  if (((ih + 1) % 16 == 0) || ih == (list->GetSize() - 1)) {
293  ctmp.Print(psname.str().c_str());
294  ctmp.Clear("D");
295  }
296  }
297 
298  psname.str("");
299  psname << Form("cdcdedx_bdcal_%s.pdf]", sfx.data());
300  ctmp.Print(psname.str().c_str());
301 }
302 
303 //------------------------------------
304 void CDCDedxBadWireAlgorithm::plotBadWireMap(const vector<double>& vbadwires, const vector<double>& vdeadwires)
305 {
306 
307  TCanvas cmap(Form("cmap_%s", m_suffix.data()), "", 800, 800);
308  cmap.SetTitle("CDC dE/dx bad wire status");
309 
310  int total = 0;
311  TH2F* hxyAll = getHistoPattern(vbadwires, "all", total);
312  hxyAll->SetTitle(Form("wire status map (%s)", m_suffix.data()));
313  setHistCosmetics(hxyAll, kGray);
314  hxyAll->Draw();
315 
316  int nbad = 0.0;
317  TH2F* hxyBad = getHistoPattern(vbadwires, "bad", nbad);
318  if (hxyBad) {
319  setHistCosmetics(hxyBad, kRed);
320  hxyBad->Draw("same");
321  }
322 
323  int ndead = 0.0;
324  TH2F* hxyDead = getHistoPattern(vdeadwires, "dead", ndead);
325  if (hxyDead) {
326  setHistCosmetics(hxyDead, kBlack);
327  hxyDead->Draw("same");
328  }
329 
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");
337  leg->Draw();
338 
339  gStyle->SetLegendTextSize(0.025);
340  TPaveText* pt = new TPaveText(-0.30, -1.47, -0.31, -1.30, "br");
341  setTextCosmetics(pt, 0.02258064);
342 
343  TText* text = pt->AddText("CDC-wire map: counter-clockwise and start from +x");
344  text->SetTextColor(kGray + 1);
345  pt->Draw("same");
346 
347  cmap.SaveAs(Form("cdcdedx_bdcal_wiremap_%s.pdf", m_suffix.data()));
348 
349  delete hxyAll;
350  delete hxyBad;
351  delete hxyDead;
352 }
353 
354 //------------------------------------
355 TH2F* CDCDedxBadWireAlgorithm::getHistoPattern(const vector<double>& inwires, const string& suffix, int& total)
356 {
357 
358  B2INFO("Creating CDCGeometryPar object");
360 
361  TH2F* temp = new TH2F(Form("temp_%s_%s", m_suffix.data(), suffix.data()), "", 2400, -1.2, 1.2, 2400, -1.2, 1.2);
362 
363  int jwire = -1;
364  total = 0;
365  for (unsigned int ilay = 0; ilay < c_maxNSenseLayers; ++ilay) {
366  for (unsigned int iwire = 0; iwire < cdcgeo.nWiresInLayer(ilay); ++iwire) {
367  jwire++;
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") {
373  total++;
374  temp->Fill(x, y);
375  } else {
376  if (count(inwires.begin(), inwires.end(), jwire)) {
377  temp->Fill(x, y);
378  total++;
379  }
380  }
381  }
382  }
383  return temp;
384 }
385 
386 //------------------------------------
387 void CDCDedxBadWireAlgorithm::plotQaPars(map<int, vector<double>>& qapars)
388 {
389 
390  string qaname[3] = {"mean", "rms", "high_fraction"};
391 
392  double linemin[3] = {m_amean* (1 - m_meanThres), m_arms* (1 - m_rmsThres), m_fracThres * 100};
393  double linemax[3] = {m_amean* (1 + m_meanThres), m_arms* (1 + m_rmsThres), m_fracThres * 100};
394 
395  for (int iqa = 0; iqa < 3; iqa++) {
396 
397  TH1D histqa(Form("%s_%s", qaname[iqa].data(), m_suffix.data()), "", c_nwireCDC, -0.5, 14335.5);
398 
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]);
402  }
403 
404  TCanvas c_pars(Form("c_pars_%d", iqa), "", 800, 600);
405  c_pars.cd();
406  gPad->SetGridy();
407 
408  histqa.SetTitle(Form("%s vs wires (%s); wire ; %s", qaname[iqa].data(), m_suffix.data(), qaname[iqa].data()));
409  histqa.SetStats(0);
410  histqa.Draw();
411 
412  TLine* lmin = new TLine(-0.5, linemin[iqa], 14335.5, linemin[iqa]);
413  lmin->SetLineColor(kRed);
414  lmin->Draw("same");
415  TLine* lmax = new TLine(-0.5, linemax[iqa], 14335.5, linemax[iqa]);
416  lmax->SetLineColor(kRed);
417  lmax->Draw("same");
418 
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()));
421 
422  delete lmax;
423  delete lmin;
424  }
425 }
426 
427 //------------------------------------
429 {
430 
431  TCanvas cstats("cstats", "cstats", 1000, 500);
432  cstats.SetBatch(kTRUE);
433  cstats.Divide(2, 1);
434 
435  cstats.cd(1);
436  auto hestats = getObjectPtr<TH1I>("hestats");
437  if (hestats) {
438  hestats->SetName(Form("htstats_%s", m_suffix.data()));
439  hestats->SetStats(0);
440  hestats->DrawCopy("");
441  }
442 
443  cstats.cd(2);
444  auto htstats = getObjectPtr<TH1I>("htstats");
445  if (htstats) {
446  htstats->SetName(Form("htstats_%s", m_suffix.data()));
447  htstats->SetStats(0);
448  htstats->DrawCopy("");
449  }
450 
451  cstats.Print(Form("cdcdedx_bdcal_qastats_%s.pdf", m_suffix.data()));
452 }
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.