Belle II Software  release-05-01-25
DQMHistAnalysisPlotOnly.cc
1 //+
2 // File : DQMHistAnalysisPlotOnly.cc
3 // Description : DQM ANalalysis, plots a list of histograms into canvases
4 //
5 // Author : B. Spruck
6 // Date : 13 - Oct - 2017
7 //-
8 
9 
10 #include <dqm/analysis/modules/DQMHistAnalysisPlotOnly.h>
11 #include <TDirectory.h>
12 #include <TROOT.h>
13 #include <TClass.h>
14 
15 using namespace std;
16 using namespace Belle2;
17 
18 //-----------------------------------------------------------------
19 // Register the Module
20 //-----------------------------------------------------------------
21 REG_MODULE(DQMHistAnalysisPlotOnly)
22 
23 //-----------------------------------------------------------------
24 // Implementation
25 //-----------------------------------------------------------------
26 
27 
30 {
31  //Parameter definition
32  addParam("HistoList", m_histlist, "histname [, canvasname]");
33  B2DEBUG(1, "DQMHistAnalysisPlotOnly: Constructor done.");
34 }
35 
36 TH1* DQMHistAnalysisPlotOnlyModule::GetHisto(TString histoname)
37 {
38  TH1* hh1;
39  hh1 = findHist(histoname.Data());
40  if (hh1 == NULL) {
41  B2DEBUG(20, "Histo " << histoname << " not in memfile");
42  // the following code sux ... is there no root function for that?
43  TDirectory* d = gROOT;
44  TString myl = histoname;
45  TString tok;
46  Ssiz_t from = 0;
47  while (myl.Tokenize(tok, from, "/")) {
48  TString dummy;
49  Ssiz_t f;
50  f = from;
51  if (myl.Tokenize(dummy, f, "/")) { // check if its the last one
52  auto e = d->GetDirectory(tok);
53  if (e) {
54  B2DEBUG(20, "Cd Dir " << tok);
55  d = e;
56  }
57  d->cd();
58  } else {
59  break;
60  }
61  }
62  TObject* obj = d->FindObject(tok);
63  if (obj != NULL) {
64  if (obj->IsA()->InheritsFrom("TH1")) {
65  B2DEBUG(20, "Histo " << histoname << " found in mem");
66  hh1 = (TH1*)obj;
67  }
68  } else {
69  B2DEBUG(20, "Histo " << histoname << " NOT found in mem");
70  }
71  }
72 
73  if (hh1 == NULL) {
74  B2DEBUG(20, "Histo " << histoname << " not found");
75  }
76 
77  return hh1;
78 }
79 
80 DQMHistAnalysisPlotOnlyModule::~DQMHistAnalysisPlotOnlyModule() { }
81 
82 void DQMHistAnalysisPlotOnlyModule::initialize()
83 {
84  B2DEBUG(20, "DQMHistAnalysisPlotOnly: initialized.");
85 
86  m_canvasList.clear();
87  for (auto& it : m_histlist) {
88  TString a = "";
89  if (it.size() == 1) {
90  a = TString(it.at(0).c_str());
91  auto p = a.Last('/');
92  if (p == kNPOS) a = "c_" + a;
93  else a.Insert(p + 1, "c_");
94  } else if (it.size() == 2) {
95  a = it.at(1).c_str();
96  } else continue;
97  TCanvas* c = new TCanvas(a);
98  m_canvasList[it.at(0)] = c;
99  }
100 }
101 
102 
103 void DQMHistAnalysisPlotOnlyModule::beginRun()
104 {
105  B2DEBUG(20, "DQMHistAnalysisPlotOnly: beginRun called.");
106 }
107 
108 void DQMHistAnalysisPlotOnlyModule::event()
109 {
110  for (auto& it : m_canvasList) {
111  gROOT->cd();
112 // TH1 *hh1 = findHist(TString(it.first.c_str()));
113  TH1* hh1 = GetHisto(it.first.c_str());
114  if (hh1 == NULL) {
115  B2DEBUG(10, "Hist " << it.first.c_str() << " not found");
116  continue;
117 // hh1 = findHistLocal(it.first);
118  }
119 
120 // if (it.second) {
121  it.second->cd();
122  if (hh1->GetDimension() == 1) {
123  hh1->Draw("hist");
124  } else if (hh1->GetDimension() == 2) {
125  hh1->Draw("colz");
126  }
127  it.second->Modified();
128  it.second->Update();
129 // }
130 
131  }
132 }
133 
134 void DQMHistAnalysisPlotOnlyModule::endRun()
135 {
136  B2DEBUG(20, "DQMHistAnalysisPlotOnly: endRun called");
137 }
138 
139 
140 void DQMHistAnalysisPlotOnlyModule::terminate()
141 {
142  B2DEBUG(20, "DQMHistAnalysisPlotOnly: terminate called");
143 }
144 
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::DQMHistAnalysisPlotOnlyModule
The module to plot a list of histograms into canvases.
Definition: DQMHistAnalysisPlotOnly.h:25
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27