Belle II Software  release-05-02-19
DQMHistAnalysis.cc
1 //+
2 // File : DQMHistAnalysisModule.cc
3 // Description : DQM histgram analysis module
4 //
5 // Author : Tomoyuki Konno, Tokyo Metropolitan Univerisity
6 // Date : 20 - Dec - 2015
7 //-
8 
9 #include <dqm/analysis/modules/DQMHistAnalysis.h>
10 #include <TROOT.h>
11 #include <TClass.h>
12 
13 
14 
15 using namespace std;
16 using namespace Belle2;
17 
18 //-----------------------------------------------------------------
19 // Register the Module
20 //-----------------------------------------------------------------
21 REG_MODULE(DQMHistAnalysis)
22 
23 //-----------------------------------------------------------------
24 // Implementation
25 //-----------------------------------------------------------------
26 
27 DQMHistAnalysisModule::ParamTypeList DQMHistAnalysisModule::g_parname;
28 DQMHistAnalysisModule::IntValueList DQMHistAnalysisModule::g_vint;
29 DQMHistAnalysisModule::FloatValueList DQMHistAnalysisModule::g_vfloat;
30 DQMHistAnalysisModule::TextList DQMHistAnalysisModule::g_text;
31 DQMHistAnalysisModule::HistList DQMHistAnalysisModule::g_hist;
32 DQMHistAnalysisModule::MonObjList DQMHistAnalysisModule::g_monObj;
33 
34 DQMHistAnalysisModule::DQMHistAnalysisModule() : Module()
35 {
36  //Set module properties
37  setDescription("Histgram Analysis module");
38 }
39 
40 
41 DQMHistAnalysisModule::~DQMHistAnalysisModule()
42 {
43 
44 }
45 
46 void DQMHistAnalysisModule::addHist(const std::string& dirname, const std::string& histname, TH1* h)
47 {
48  if (dirname.size() > 0) {
49  g_hist.insert(HistList::value_type(dirname + "/" + histname, h));
50  } else {
51  g_hist.insert(HistList::value_type(histname, h));
52  }
53 }
54 
56 {
57  if (g_monObj.find(objName) != g_monObj.end()) {
58  if (g_monObj[objName]) {
59  return g_monObj[objName];
60  } else {
61  B2WARNING("MonitoringObject " << objName << " listed as being in memfile but points to nowhere. New Object will be made.");
62  g_monObj.erase(objName);
63  }
64  }
65 
66  MonitoringObject* obj = new MonitoringObject(objName);
67  g_monObj.insert(MonObjList::value_type(objName, obj));
68  return obj;
69 }
70 
71 
73 {
74  return g_hist;
75 }
76 
78 {
79  return g_monObj;
80 }
81 
82 
83 TH1* DQMHistAnalysisModule::findHist(const std::string& histname)
84 {
85  if (g_hist.find(histname) != g_hist.end()) {
86  if (g_hist[histname]) {
87  //Want to search elsewhere if null-pointer saved in map
88  return g_hist[histname];
89  } else {
90  B2ERROR("Histogram " << histname << " listed as being in memfile but points to nowhere.");
91  }
92  }
93  B2INFO("Histogram " << histname << " not in memfile.");
94 
95  //Histogram not in list, search in memory for it
96  gROOT->cd();
97 
98  //Following the path to the histogram
99  TDirectory* d = gROOT;
100  TString myl = histname;
101  TString tok;
102  Ssiz_t from = 0;
103  while (myl.Tokenize(tok, from, "/")) {
104  TString dummy;
105  Ssiz_t f;
106  f = from;
107  if (myl.Tokenize(dummy, f, "/")) { // check if its the last one
108  auto e = d->GetDirectory(tok);
109  if (e) {
110  B2INFO("Cd Dir " << tok);
111  d = e;
112  }
113  d->cd();
114  } else {
115  break;
116  }
117  }
118 
119  // This code assumes that the histograms address does NOT change between initialization and any later event
120  // This assumption seems to be reasonable for TFiles and in-memory objects
121  // BUT this means => Analysis moules MUST NEVER create a histogram with already existing name NOR delete any histogram
122  TH1* found_hist = findHist(d, tok);
123  if (found_hist) {
124  g_hist[histname] = found_hist;//Can't use addHist as we want to overwrite invalid entries
125  }
126  return found_hist;
127 
128 }
129 
130 TH1* DQMHistAnalysisModule::findHist(const std::string& dirname, const std::string& histname)
131 {
132  if (dirname.size() > 0) {
133  return findHist(dirname + "/" + histname);
134  }
135  return findHist(histname);
136 }
137 
138 
139 TH1* DQMHistAnalysisModule::findHist(const TDirectory* histdir, const TString& histname)
140 {
141  TObject* obj = histdir->FindObject(histname);
142  if (obj != NULL) {
143  if (obj->IsA()->InheritsFrom("TH1")) {
144  B2INFO("Histogram " << histname << " found in mem");
145  return (TH1*)obj;
146  }
147  } else {
148  B2INFO("Histogram " << histname << " NOT found in mem");
149  }
150  return NULL;
151 }
152 
154 {
155  if (g_monObj.find(objName) != g_monObj.end()) {
156  if (g_monObj[objName]) {
157  //Want to search elsewhere if null-pointer saved in map
158  return g_monObj[objName];
159  } else {
160  B2ERROR("MonitoringObject " << objName << " listed as being in memfile but points to nowhere.");
161  }
162  }
163  B2INFO("MonitoringObject " << objName << " not in memfile.");
164  return NULL;
165 }
166 
167 
168 
169 void DQMHistAnalysisModule::setIntValue(const std::string& parname, int vint)
170 {
171  if (g_parname.find(parname) == g_parname.end() && g_vint.find(parname) == g_vint.end()) {
172  g_parname.insert(ParamTypeList::value_type(parname, c_ParamINT));
173  g_vint.insert(IntValueList::value_type(parname, vint));
174  } else if (g_vint.find(parname) == g_vint.end()) {
175  B2ERROR(parname + " is already registered as non-int data type");
176  } else {
177  g_vint[parname] = vint;
178  }
179 }
180 
181 void DQMHistAnalysisModule::setFloatValue(const std::string& parname, float vfloat)
182 {
183  if (g_parname.find(parname) == g_parname.end() && g_vfloat.find(parname) == g_vfloat.end()) {
184  g_parname.insert(ParamTypeList::value_type(parname, c_ParamFLOAT));
185  g_vfloat.insert(FloatValueList::value_type(parname, vfloat));
186  } else if (g_vfloat.find(parname) == g_vfloat.end()) {
187  B2ERROR(parname + " is already registered as non-float data type");
188  } else {
189  g_vfloat[parname] = vfloat;
190  }
191 }
192 
193 void DQMHistAnalysisModule::setText(const std::string& parname, const std::string& text)
194 {
195  if (g_parname.find(parname) == g_parname.end() && g_text.find(parname) == g_text.end()) {
196  g_parname.insert(ParamTypeList::value_type(parname, c_ParamTEXT));
197  g_text.insert(TextList::value_type(parname, text));
198  } else if (g_text.find(parname) == g_text.end()) {
199  B2ERROR(parname + " is already registered as non-text data type");
200  } else {
201  g_text[parname] = text;
202  }
203 }
204 
Belle2::DQMHistAnalysisModule::getMonObjList
static const MonObjList & getMonObjList()
Get the list of MonitoringObjects.
Definition: DQMHistAnalysis.cc:77
Belle2::DQMHistAnalysisModule::getHistList
static const HistList & getHistList()
Get the list of the histograms.
Definition: DQMHistAnalysis.cc:72
Belle2::Module::setDescription
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:216
Belle2::DQMHistAnalysisModule::setFloatValue
static void setFloatValue(const std::string &parname, float vfloat)
Set the float value of the parameter.
Definition: DQMHistAnalysis.cc:181
Belle2::DQMHistAnalysisModule::g_parname
static ParamTypeList g_parname
The list of module parameter types.
Definition: DQMHistAnalysis.h:78
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::DQMHistAnalysisModule::setIntValue
static void setIntValue(const std::string &parname, int vint)
Set the integer value of the parameter.
Definition: DQMHistAnalysis.cc:169
Belle2::DQMHistAnalysisModule::findHist
static TH1 * findHist(const std::string &histname)
Find histogram.
Definition: DQMHistAnalysis.cc:83
Belle2::DQMHistAnalysisModule::c_ParamTEXT
@ c_ParamTEXT
The string type for module parameter.
Definition: DQMHistAnalysis.h:45
Belle2::DQMHistAnalysisModule::TextList
std::map< std::string, std::string > TextList
The type of list of string module parameter.
Definition: DQMHistAnalysis.h:62
Belle2::DQMHistAnalysisModule::g_hist
static HistList g_hist
The list of histograms.
Definition: DQMHistAnalysis.h:94
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::DQMHistAnalysisModule::IntValueList
std::map< std::string, int > IntValueList
The type of list of integer module parameter.
Definition: DQMHistAnalysis.h:54
Belle2::DQMHistAnalysisModule::g_vfloat
static FloatValueList g_vfloat
The list of float module parameter.
Definition: DQMHistAnalysis.h:86
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisModule::HistList
std::map< std::string, TH1 * > HistList
The type of list of histograms.
Definition: DQMHistAnalysis.h:66
Belle2::DQMHistAnalysisModule::g_monObj
static MonObjList g_monObj
The list of MonitoringObjects.
Definition: DQMHistAnalysis.h:98
Belle2::DQMHistAnalysisModule::FloatValueList
std::map< std::string, float > FloatValueList
The type of list of float module parameter.
Definition: DQMHistAnalysis.h:58
Belle2::DQMHistAnalysisModule::setText
static void setText(const std::string &parname, const std::string &text)
Set the string value of the parameter.
Definition: DQMHistAnalysis.cc:193
Belle2::DQMHistAnalysisModule::g_text
static TextList g_text
The list of string module parameter.
Definition: DQMHistAnalysis.h:90
Belle2::DQMHistAnalysisModule::c_ParamFLOAT
@ c_ParamFLOAT
The float type for module parameter.
Definition: DQMHistAnalysis.h:41
Belle2::DQMHistAnalysisModule::g_vint
static IntValueList g_vint
The list of integer module parameter.
Definition: DQMHistAnalysis.h:82
Belle2::DQMHistAnalysisModule::findMonitoringObject
static MonitoringObject * findMonitoringObject(const std::string &objName)
Find MonitoringObject.
Definition: DQMHistAnalysis.cc:153
Belle2::DQMHistAnalysisModule::ParamTypeList
std::map< std::string, EParamType > ParamTypeList
The type of list of module parameter types.
Definition: DQMHistAnalysis.h:50
Belle2::DQMHistAnalysisModule::MonObjList
std::map< std::string, MonitoringObject * > MonObjList
The type of list of MonitoringObjects.
Definition: DQMHistAnalysis.h:70
Belle2::DQMHistAnalysisModule::getMonitoringObject
static MonitoringObject * getMonitoringObject(const std::string &histname)
Get MonitoringObject with given name (new object is created if non-existing)
Definition: DQMHistAnalysis.cc:55
Belle2::MonitoringObject
MonitoringObject is a basic object to hold data for the run-dependency monitoring Run summary TCanvas...
Definition: MonitoringObject.h:41
Belle2::DQMHistAnalysisModule::addHist
static void addHist(const std::string &dirname, const std::string &histname, TH1 *h)
Add histogram.
Definition: DQMHistAnalysis.cc:46
Belle2::DQMHistAnalysisModule::c_ParamINT
@ c_ParamINT
The integer type for module parameter.
Definition: DQMHistAnalysis.h:37