Belle II Software  release-06-00-14
DQMHistAnalysisPXDTrackCharge.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 // File : DQMHistAnalysisPXDTrackCharge.cc
10 // Description : Analysis of PXD Cluster Charge
11 //-
12 
13 
14 #include <dqm/analysis/modules/DQMHistAnalysisPXDTrackCharge.h>
15 #include <TROOT.h>
16 #include <TStyle.h>
17 #include <TClass.h>
18 #include <TLatex.h>
19 #include <vxd/geometry/GeoCache.h>
20 #include <TKey.h>
21 
22 #include <RooDataHist.h>
23 #include <RooAbsPdf.h>
24 #include <RooPlot.h>
25 
26 
27 using namespace std;
28 using namespace Belle2;
29 
30 //-----------------------------------------------------------------
31 // Register the Module
32 //-----------------------------------------------------------------
33 REG_MODULE(DQMHistAnalysisPXDTrackCharge)
34 
35 //-----------------------------------------------------------------
36 // Implementation
37 //-----------------------------------------------------------------
38 
41 {
42  // This module CAN NOT be run in parallel!
43 
44  //Parameter definition
45  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of Histogram dir", std::string("PXDER"));
46  addParam("RangeLow", m_rangeLow, "Lower border for fit", 20.);
47  addParam("RangeHigh", m_rangeHigh, "High border for fit", 80.);
48 // addParam("PeakBefore", m_peakBefore, "Range for fit before peak (positive)", 5.);
49 // addParam("PeakAfter", m_peakAfter, "Range for after peak", 40.);
50  addParam("PVPrefix", m_pvPrefix, "PV Prefix", std::string("DQM:PXD:TrackCharge:"));
51  addParam("useEpics", m_useEpics, "useEpics", true);
52  addParam("RefHistoFile", m_refFileName, "Reference histrogram file name", std::string("refHisto.root"));
53  addParam("ColorAlert", m_color, "Whether to show the color alert", true);
54 
55  B2DEBUG(99, "DQMHistAnalysisPXDTrackCharge: Constructor done.");
56 }
57 
58 DQMHistAnalysisPXDTrackChargeModule::~DQMHistAnalysisPXDTrackChargeModule()
59 {
60 #ifdef _BELLE2_EPICS
61  if (m_useEpics) {
62  if (ca_current_context()) ca_context_destroy();
63  }
64 #endif
65 }
66 
67 void DQMHistAnalysisPXDTrackChargeModule::initialize()
68 {
69  B2DEBUG(99, "DQMHistAnalysisPXDTrackCharge: initialized.");
70 
71  m_refFile = NULL;
72  if (m_refFileName != "") {
73  m_refFile = new TFile(m_refFileName.data());
74  }
75 
76  m_monObj = getMonitoringObject("pxd");
77  const VXD::GeoCache& geo = VXD::GeoCache::getInstance();
78 
79  m_rfws = new RooWorkspace("w");
80  m_rfws->factory("Landau::landau(x[0,100],ml[20,10,50],sl[5,1,30])");
81  m_rfws->factory("Gaussian::gauss(x,mg[0],sg[2,0.1,10])");
82  m_rfws->factory("FCONV::lxg(x,landau,gauss)");
83 
84  m_x = m_rfws->var("x");
85  m_x->setRange("signal", m_rangeLow, m_rangeHigh);
86 
87  // collect the list of all PXD Modules in the geometry here
88  std::vector<VxdID> sensors = geo.getListOfSensors();
89  for (VxdID& aVxdID : sensors) {
90  VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
91  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
92  m_PXDModules.push_back(aVxdID); // reorder, sort would be better
93 
94  std::string name = "PXD_Track_Cluster_Charge_" + (std::string)aVxdID;
95  std::replace(name.begin(), name.end(), '.', '_');
96  m_cChargeMod[aVxdID] = new TCanvas((m_histogramDirectoryName + "/c_Fit_" + name).data());
97  if (aVxdID == VxdID("1.5.1")) {
98  for (int s = 0; s < 6; s++) {
99  for (int d = 0; d < 4; d++) {
100  m_cChargeModASIC[aVxdID][s][d] = new TCanvas((m_histogramDirectoryName + "/c_Fit_" + name + Form("_s%d_d%d", s + 1, d + 1)).data());
101  }
102  }
103  m_hChargeModASIC2d[aVxdID] = new TH2F(("hPXD_TCChargeMPV_" + name).data(),
104  ("PXD TCCharge MPV " + name + ";Switcher;DCD;MPV").data(),
105  6, 0.5, 6.5, 4, 0.5, 4.5);
106  m_cChargeModASIC2d[aVxdID] = new TCanvas((m_histogramDirectoryName + "/c_TCCharge_MPV_" + name).data());
107  }
108  }
109  std::sort(m_PXDModules.begin(), m_PXDModules.end()); // back to natural order
110 
111  gROOT->cd(); // this seems to be important, or strange things happen
112 
113  m_cTrackedClusters = new TCanvas((m_histogramDirectoryName + "/c_TrackedClusters").data());
114  m_hTrackedClusters = new TH1F("hPXDTrackedClusters", "PXD Tracked Clusters/Event;Module", 40, 0, 40);
115  m_hTrackedClusters->Draw();
116  auto ax = m_hTrackedClusters->GetXaxis();
117  if (ax) {
118  ax->Set(m_PXDModules.size(), 0, m_PXDModules.size());
119  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
120  TString ModuleName = (std::string)m_PXDModules[i];
121  ax->SetBinLabel(i + 1, ModuleName);
122  }
123  } else B2ERROR("no axis");
124 
125  m_cCharge = new TCanvas((m_histogramDirectoryName + "/c_TrackCharge").data());
126  m_monObj->addCanvas(m_cCharge);
127 
128  m_gCharge = new TGraphErrors();
129  m_gCharge->SetName("Track_Cluster_Charge");
130  m_gCharge->SetTitle("Track Cluster Charge");
131 
133  m_line_up = new TLine(0, 10, m_PXDModules.size(), 10);
134  m_line_mean = new TLine(0, 16, m_PXDModules.size(), 16);
135  m_line_low = new TLine(0, 3, m_PXDModules.size(), 3);
136  m_line_up->SetHorizontal(true);
137  m_line_up->SetLineColor(kMagenta);// Green
138  m_line_up->SetLineWidth(3);
139  m_line_up->SetLineStyle(7);
140  m_line_mean->SetHorizontal(true);
141  m_line_mean->SetLineColor(kGreen);// Black
142  m_line_mean->SetLineWidth(3);
143  m_line_mean->SetLineStyle(4);
144  m_line_low->SetHorizontal(true);
145  m_line_low->SetLineColor(kMagenta);
146  m_line_low->SetLineWidth(3);
147  m_line_low->SetLineStyle(7);
148 
149  m_fMean = new TF1("f_Mean", "pol0", 0, m_PXDModules.size());
150  m_fMean->SetParameter(0, 50);
151  m_fMean->SetLineColor(kYellow);
152  m_fMean->SetLineWidth(3);
153  m_fMean->SetLineStyle(7);
154  m_fMean->SetNpx(m_PXDModules.size());
155  m_fMean->SetNumberFitPoints(m_PXDModules.size());
156 
157 #ifdef _BELLE2_EPICS
158  if (m_useEpics) {
159  if (!ca_current_context()) SEVCHK(ca_context_create(ca_disable_preemptive_callback), "ca_context_create");
160  mychid.resize(3);
161  SEVCHK(ca_create_channel((m_pvPrefix + "Mean").data(), NULL, NULL, 10, &mychid[0]), "ca_create_channel failure");
162  SEVCHK(ca_create_channel((m_pvPrefix + "Diff").data(), NULL, NULL, 10, &mychid[1]), "ca_create_channel failure");
163  SEVCHK(ca_create_channel((m_pvPrefix + "Status").data(), NULL, NULL, 10, &mychid[2]), "ca_create_channel failure");
164  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
165  }
166 #endif
167 }
168 
169 
170 void DQMHistAnalysisPXDTrackChargeModule::beginRun()
171 {
172  B2DEBUG(99, "DQMHistAnalysisPXDTrackCharge: beginRun called.");
173 
174  m_cCharge->Clear();
175 }
176 
177 void DQMHistAnalysisPXDTrackChargeModule::event()
178 {
179  if (!m_cCharge) return;
180 
181  gStyle->SetOptStat(0);
182  gStyle->SetStatStyle(1);
183  gStyle->SetOptDate(22);// Date and Time in Bottom Right, does no work
184 
185  bool enough = false;
186 
187 // auto landau = m_rfws->pdf("landau");
188 // auto gauss = m_rfws->pdf("gauss");
189  auto model = m_rfws->pdf("lxg");
190 
191  auto ml = m_rfws->var("ml");
192 // auto sl = m_rfws->var("sl");
193 // auto mg = m_rfws->var("mg");
194 // auto sg = m_rfws->var("sg");
195 
196  {
197  m_cTrackedClusters->Clear();
198  m_cTrackedClusters->cd();
199  m_hTrackedClusters->Reset();
200 
201  std::string name = "Tracked_Clusters"; // new name
202  TH1* hh2 = findHist(m_histogramDirectoryName, "PXD_Tracked_Clusters");
203  if (hh2) {
204  auto scale = hh2->GetBinContent(0);// overflow misused as event counter!
205  if (scale > 0) {
206  int j = 1;
207  for (int i = 0; i < 64; i++) {
208  auto layer = (((i >> 5) & 0x1) + 1);
209  auto ladder = ((i >> 1) & 0xF);
210  auto sensor = ((i & 0x1) + 1);
211 
212  auto id = Belle2::VxdID(layer, ladder, sensor);
213  // Check if sensor exist
214  if (Belle2::VXD::GeoCache::getInstance().validSensorID(id)) {
215  m_hTrackedClusters->SetBinContent(j, hh2->GetBinContent(i + 1) / scale);
216  j++;
217  }
218  }
219  }
220  m_hTrackedClusters->SetName(name.data());
221  m_hTrackedClusters->SetTitle("Tracked Clusters/Event");
222  m_hTrackedClusters->SetFillColor(kWhite);
223  m_hTrackedClusters->SetStats(kFALSE);
224  m_hTrackedClusters->SetLineStyle(1);// 2 or 3
225  m_hTrackedClusters->SetLineColor(kBlue);
226  m_hTrackedClusters->Draw("hist");
227 
228  TH1* href2 = GetHisto("ref/" + m_histogramDirectoryName + "/" + name);
229 
230  if (href2) {
231  href2->SetLineStyle(3);// 2 or 3
232  href2->SetLineColor(kBlack);
233  href2->Draw("same,hist");
234  }
235 
236  auto tt = new TLatex(5.5, 0, " 1.3.2 Module is excluded, please ignore");
237  tt->SetTextAngle(90);// Rotated
238  tt->SetTextAlign(12);// Centered
239  tt->Draw();
240  }
241  }
242 
243  m_gCharge->Set(0);
244 
245  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
246  TCanvas* canvas = m_cChargeMod[m_PXDModules[i]];
247  if (canvas == nullptr) continue;
248 
249  std::string name = "PXD_Track_Cluster_Charge_" + (std::string)m_PXDModules[i];
250  std::replace(name.begin(), name.end(), '.', '_');
251 
252  canvas->cd();
253  canvas->Clear();
254 
255  TH1* hh1 = findHist(m_histogramDirectoryName, name);
256  if (hh1) {
257 
258  if (hh1->GetEntries() > 50) {
259 
260  auto hdata = new RooDataHist(hh1->GetName(), hh1->GetTitle(), *m_x, (const TH1*) hh1);
261  auto plot = m_x->frame(RooFit::Title(hh1->GetTitle()));
262  /*auto r =*/ model->fitTo(*hdata, RooFit::Range("signal"));
263 
264  model->paramOn(plot, RooFit::Format("NELU", RooFit::AutoPrecision(2)), RooFit::Layout(0.6, 0.9, 0.9));
265  hdata->plotOn(plot, RooFit::LineColor(kBlue)/*, RooFit::Range("plot"), RooFit::NormRange("signal")*/);
266  model->plotOn(plot, RooFit::LineColor(kRed), RooFit::Range("signal"), RooFit::NormRange("signal"));
267 
268  plot->Draw("");
269 
270 // model->Print("");
271 // ml->Print("");
272 // sl->Print("");
273 // mg->Print("");
274 // sg->Print("");
275 // cout << "ZZZ , " << Form("%d%02d%d ,", std::get<0>(t), std::get<1>(t), std::get<2>(t)) << ml->getValV() << "," << ml->getError() << "," << sl->getValV() << "," << sl->getError() << "," << sg->getValV() << "," << sg->getError() << endl;
276 
277 
278  int p = m_gCharge->GetN();
279  m_gCharge->SetPoint(p, i + 0.49, ml->getValV());
280  m_gCharge->SetPointError(p, 0.1, ml->getError()); // error in x is useless
281  m_monObj->setVariable(("trackcharge_" + (std::string)m_PXDModules[i]).c_str(), ml->getValV(), ml->getError());
282  }
283 
284  TH1* hist2 = GetHisto("ref/" + m_histogramDirectoryName + "/" + name);
285 
286  if (hist2) {
287 // B2INFO("Draw Normalized " << hist2->GetName());
288  hist2->SetLineStyle(3);// 2 or 3
289  hist2->SetLineColor(kBlack);
290 
291 // TIter nextkey(canvas->GetListOfPrimitives());
292 // TObject* obj = NULL;
293 // while ((obj = (TObject*)nextkey())) {
294 // if (obj->IsA()->InheritsFrom("TH1")) {
295 // if (string(obj->GetName()) == string(hist2->GetName())) {
296 // delete obj;
297 // }
298 // }
299 // }
300 
301  canvas->cd();
302 
303  // if draw normalized
304  TH1* h = (TH1*)hist2->Clone(); // Annoying ... Maybe an memory leak? TODO
305  if (abs(hist2->GetEntries()) > 0) h->Scale(hh1->GetEntries() / hist2->GetEntries());
306 
307  h->SetStats(kFALSE);
308  h->Draw("same,hist");
309  }
310 
311  // add coloring, cuts? based on fit, compare with ref?
312  canvas->Pad()->SetFrameFillColor(10);
313  if (m_color) {
314  if (hh1->GetEntries() < 100) {
315  // not enough Entries
316  canvas->Pad()->SetFillColor(kGray);
317  } else {
318  canvas->Pad()->SetFillColor(kGreen);
319  }
320  } else {
321  canvas->Pad()->SetFillColor(kWhite);// White
322  }
323 
324  canvas->Modified();
325  canvas->Update();
326 
327  // means if ANY plot is > 100 entries, all plots are assumed to be o.k.
328  if (hh1->GetEntries() >= 1000) enough = true;
329  }
330  }
331 
332  // now loop per module over asics pairs (1.5.1)
333  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
334 // TCanvas* canvas = m_cChargeMod[m_PXDModules[i]];
335  VxdID& aVxdID = m_PXDModules[i];
336 
337  if (m_hChargeModASIC2d[aVxdID]) m_hChargeModASIC2d[aVxdID]->Reset();
338  if (m_cChargeModASIC2d[aVxdID]) m_cChargeModASIC2d[aVxdID]->Clear();
339 
340  for (int s = 1; s <= 6; s++) {
341  for (int d = 1; d <= 4; d++) {
342  std::string name = "PXD_Track_Cluster_Charge_" + (std::string)m_PXDModules[i] + Form("_sw%d_dcd%d", s, d);
343  std::replace(name.begin(), name.end(), '.', '_');
344 
345  if (m_cChargeModASIC[aVxdID][s - 1][d - 1]) {
346  m_cChargeModASIC[aVxdID][s - 1][d - 1]->Clear();
347  m_cChargeModASIC[aVxdID][s - 1][d - 1]->cd();
348  }
349 
350  TH1* hh1 = findHist(m_histogramDirectoryName, name);
351  if (hh1) {
352  double mpv = 0.0;
353  if (hh1->GetEntries() > 50) {
354  auto hdata = new RooDataHist(hh1->GetName(), hh1->GetTitle(), *m_x, (const TH1*) hh1);
355  auto plot = m_x->frame(RooFit::Title(hh1->GetTitle()));
356  /*auto r =*/ model->fitTo(*hdata, RooFit::Range("signal"));
357 
358  if (m_cChargeModASIC[aVxdID][s - 1][d - 1]) {
359  model->paramOn(plot, RooFit::Format("NELU", RooFit::AutoPrecision(2)), RooFit::Layout(0.6, 0.9, 0.9));
360  hdata->plotOn(plot, RooFit::LineColor(kBlue)/*, RooFit::Range("plot"), RooFit::NormRange("signal")*/);
361  model->plotOn(plot, RooFit::LineColor(kRed), RooFit::Range("signal"), RooFit::NormRange("signal"));
362  }
363  plot->Draw("");
364 
365  mpv = ml->getValV();
366  }
367 
368  if (m_hChargeModASIC2d[aVxdID]) {
369  if (mpv > 0.0) m_hChargeModASIC2d[aVxdID]->Fill(s, d, mpv); // TODO check what is s, d
370  }
371  }
372  }
373  }
374 
375  // Overview map of ASCI combinations
376  if (m_hChargeModASIC2d[aVxdID] && m_cChargeModASIC2d[aVxdID]) {
377  m_cChargeModASIC2d[aVxdID]->cd();
378  m_hChargeModASIC2d[aVxdID]->Draw("colz");
379  }
380  }
381 
382  m_cCharge->cd();
383  m_cCharge->Clear();
384  m_gCharge->SetMinimum(0);
385  m_gCharge->SetMaximum(70);
386  auto ax = m_gCharge->GetXaxis();
387  if (ax) {
388  ax->Set(m_PXDModules.size(), 0, m_PXDModules.size());
389  for (unsigned int i = 0; i < m_PXDModules.size(); i++) {
390  TString ModuleName = (std::string)m_PXDModules[i];
391  ax->SetBinLabel(i + 1, ModuleName);
392  }
393  } else B2ERROR("no axis");
394 
395  m_gCharge->SetLineColor(4);
396  m_gCharge->SetLineWidth(2);
397  m_gCharge->SetMarkerStyle(8);
398  m_gCharge->Draw("AP");
399  {
400  auto tt = new TLatex(5.5, 0, " 1.3.2 Module is excluded, please ignore");
401  tt->SetTextAngle(90);// Rotated
402  tt->SetTextAlign(12);// Centered
403  tt->Draw();
404  }
405  m_cCharge->cd(0);
406  m_cCharge->Modified();
407  m_cCharge->Update();
408 
409  double data = 0;
410  double diff = 0;
411  if (m_gCharge && enough) {
412 // double currentMin, currentMax;
413  m_gCharge->Fit(m_fMean, "R");
414  double mean = m_gCharge->GetMean(2);
415  double maxi = mean + 15;
416  double mini = mean - 15;
417  m_line_up->SetY1(maxi);
418  m_line_up->SetY2(maxi);
419  m_line_mean->SetY1(mean);
420  m_line_mean->SetY2(mean);
421  m_line_low->SetY1(mini);
422  m_line_low->SetY2(mini);
423  data = mean; // m_fMean->GetParameter(0); // we are more interessted in the maximum deviation from mean
424  // m_gCharge->GetMinimumAndMaximum(currentMin, currentMax);
425  diff = m_gCharge->GetRMS(2);// RMS of Y
426  // better, max deviation as fabs(data - currentMin) > fabs(currentMax - data) ? fabs(data - currentMin) : fabs(currentMax - data);
427  m_line_up->Draw();
428  m_line_mean->Draw();
429  m_line_low->Draw();
430 
431  m_monObj->setVariable("trackcharge", mean, diff);
432  }
433 
434 #ifdef _BELLE2_EPICS
435  if (m_useEpics) {
436  SEVCHK(ca_put(DBR_DOUBLE, mychid[0], (void*)&data), "ca_set failure");
437  SEVCHK(ca_put(DBR_DOUBLE, mychid[1], (void*)&diff), "ca_set failure");
438  }
439 #endif
440 
441  int status = 0;
442 
443  if (!enough) {
444  // not enough Entries
445  m_cCharge->Pad()->SetFillColor(kGray);// Magenta or Gray
446  status = 0; // default
447  } else {
449  if (fabs(data - 30.) > 20. || diff > 12) {
450  m_cCharge->Pad()->SetFillColor(kRed);// Red
451  status = 4;
452  } else if (fabs(data - 30) > 15. || diff > 8) {
453  m_cCharge->Pad()->SetFillColor(kYellow);// Yellow
454  status = 3;
455  } else {
456  m_cCharge->Pad()->SetFillColor(kGreen);// Green
457  status = 2;
458  }
459 
460  // FIXME overwrite for now
461  // m_cCharge->Pad()->SetFillColor(kGreen);// Green
462 
463  }
464 
465 #ifdef _BELLE2_EPICS
466  if (m_useEpics) {
467  SEVCHK(ca_put(DBR_INT, mychid[2], (void*)&status), "ca_set failure");
468  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
469  }
470 #endif
471 
472 }
473 
474 void DQMHistAnalysisPXDTrackChargeModule::endRun()
475 {
476  B2DEBUG(99, "DQMHistAnalysisPXDTrackCharge : endRun called");
477 }
478 
479 
480 void DQMHistAnalysisPXDTrackChargeModule::terminate()
481 {
482  B2DEBUG(99, "DQMHistAnalysisPXDTrackCharge: terminate called");
483  // should delete canvas here, maybe hist, too? Who owns it?
484 #ifdef _BELLE2_EPICS
485  if (m_useEpics) {
486  for (auto m : mychid) SEVCHK(ca_clear_channel(m), "ca_clear_channel failure");
487  SEVCHK(ca_pend_io(5.0), "ca_pend_io failure");
488  }
489 #endif
490  if (m_refFile) delete m_refFile;
491 
492 }
493 
494 TH1* DQMHistAnalysisPXDTrackChargeModule::GetHisto(TString histoname)
495 {
496  TH1* hh1 = nullptr;
497  gROOT->cd();
498 // hh1 = findHist(histoname.Data());
499  // cppcheck-suppress knownConditionTrueFalse
500  if (hh1 == NULL) {
501  B2DEBUG(20, "findHisto failed " << histoname << " not in memfile");
502 
503  // first search reference root file ... if ther is one
504  if (m_refFile && m_refFile->IsOpen()) {
505  TDirectory* d = m_refFile;
506  TString myl = histoname;
507  TString tok;
508  Ssiz_t from = 0;
509  B2DEBUG(20, myl);
510  while (myl.Tokenize(tok, from, "/")) {
511  TString dummy;
512  Ssiz_t f;
513  f = from;
514  if (myl.Tokenize(dummy, f, "/")) { // check if its the last one
515  auto e = d->GetDirectory(tok);
516  if (e) {
517  B2DEBUG(20, "Cd Dir " << tok << " from " << d->GetPath());
518  d = e;
519  } else {
520  B2DEBUG(20, "cd failed " << tok << " from " << d->GetPath());
521  }
522  } else {
523  break;
524  }
525  }
526  TObject* obj = d->FindObject(tok);
527  if (obj != NULL) {
528  if (obj->IsA()->InheritsFrom("TH1")) {
529  B2DEBUG(20, "Histo " << histoname << " found in ref file");
530  hh1 = (TH1*)obj;
531  } else {
532  B2DEBUG(20, "Histo " << histoname << " found in ref file but wrong type");
533  }
534  } else {
535  // seems find will only find objects, not keys, thus get the object on first access
536  TIter next(d->GetListOfKeys());
537  TKey* key;
538  while ((key = (TKey*)next())) {
539  TObject* obj2 = key->ReadObj() ;
540  if (obj2->InheritsFrom("TH1")) {
541  if (obj2->GetName() == tok) {
542  hh1 = (TH1*)obj2;
543  B2DEBUG(20, "Histo " << histoname << " found as key -> readobj");
544  break;
545  }
546  }
547  }
548  if (hh1 == NULL) B2DEBUG(20, "Histo " << histoname << " NOT found in ref file " << tok);
549  }
550  }
551 
552  if (hh1 == NULL) {
553  B2DEBUG(20, "Histo " << histoname << " not in memfile or ref file");
554 
555  TDirectory* d = gROOT;
556  TString myl = histoname;
557  TString tok;
558  Ssiz_t from = 0;
559  while (myl.Tokenize(tok, from, "/")) {
560  TString dummy;
561  Ssiz_t f;
562  f = from;
563  if (myl.Tokenize(dummy, f, "/")) { // check if its the last one
564  auto e = d->GetDirectory(tok);
565  if (e) {
566  B2DEBUG(20, "Cd Dir " << tok);
567  d = e;
568  } else B2DEBUG(20, "cd failed " << tok);
569  d->cd();
570  } else {
571  break;
572  }
573  }
574  TObject* obj = d->FindObject(tok);
575  if (obj != NULL) {
576  if (obj->IsA()->InheritsFrom("TH1")) {
577  B2DEBUG(20, "Histo " << histoname << " found in mem");
578  hh1 = (TH1*)obj;
579  }
580  } else {
581  B2DEBUG(20, "Histo " << histoname << " NOT found in mem");
582  }
583  }
584  }
585 
586  if (hh1 == NULL) {
587  B2DEBUG(20, "Histo " << histoname << " not found");
588  }
589 
590  return hh1;
591 }
The base class for the histogram analysis module.
DQM Histogram Analysis for PXD Cluster Charge.
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition: GeoCache.cc:58
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:66
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:213
Base class to provide Sensor Information for PXD and SVD.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.