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