Belle II Software  release-08-01-10
DQMHistAnalysisInputRootFile.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 : DQMHistAnalysisInputRootFile.cc
10 // Description : Module for offline testing of histogram analysis code.
11 // Root file containing DQM histograms can be used as input.
12 //-
13 
14 
15 #include <dqm/analysis/modules/DQMHistAnalysisInputRootFile.h>
16 
17 #include <TROOT.h>
18 #include <TKey.h>
19 
20 #include <boost/regex.hpp>
21 #include <boost/algorithm/string/replace.hpp>
22 #include <iostream>
23 
24 using namespace Belle2;
25 
26 //-----------------------------------------------------------------
27 // Register the Module
28 //-----------------------------------------------------------------
29 REG_MODULE(DQMHistAnalysisInputRootFile);
30 
31 //-----------------------------------------------------------------
32 // Implementation
33 //-----------------------------------------------------------------
34 
37 {
38  //Parameter definition
39  addParam("FileList", m_fileList, "List of input files", std::vector<std::string> {"input_histo.root"});
40  addParam("SelectHistograms", m_histograms, "List of histogram name patterns, empty for all. Support wildcard matching (* and ?).",
41  std::vector<std::string>());
42  addParam("Experiment", m_expno, "Experiment Nr", 7u);
43  addParam("RunList", m_runList, "Run Number List", std::vector<unsigned int> {1u});
44  addParam("EventsList", m_eventsList, "Number of events for each run", std::vector<unsigned int> {10u});
45  addParam("RunType", m_runType, "Run Type override", std::string(""));
46  addParam("EventFilled", m_fillEvent, "Event override", 0);
47  addParam("EventInterval", m_interval, "Time between events (seconds)", 20u);
48  addParam("NullHistogramMode", m_nullHistoMode, "Test mode for null histograms", false);
49  B2DEBUG(1, "DQMHistAnalysisInputRootFile: Constructor done.");
50 }
51 
53 {
54  if (m_file != nullptr) delete m_file;
55  if (m_fileList.size() == 0) B2ERROR("File list is empty.");
56  if (m_runList.size() == 0) B2ERROR("Run list is empty.");
57  if (m_eventsList.size() == 0) B2ERROR("Events list is empty.");
58  if (m_runList.size() != m_eventsList.size()) B2ERROR("Run list does not have the same size as events list.");
59  if (m_runList.size() != m_fileList.size()) B2ERROR("Run list does not have the same size as file list.");
60  m_run_idx = 0;
61  m_file = new TFile(m_fileList[m_run_idx].c_str());
62  m_eventMetaDataPtr.registerInDataStore();
63  B2INFO("DQMHistAnalysisInputRootFile: initialized.");
64 }
65 
66 bool DQMHistAnalysisInputRootFileModule::hnamePatternMatch(std::string pattern, std::string text)
67 {
68  boost::replace_all(pattern, "\\", "\\\\");
69  boost::replace_all(pattern, "^", "\\^");
70  boost::replace_all(pattern, ".", "\\.");
71  boost::replace_all(pattern, "$", "\\$");
72  boost::replace_all(pattern, "|", "\\|");
73  boost::replace_all(pattern, "(", "\\(");
74  boost::replace_all(pattern, ")", "\\)");
75  boost::replace_all(pattern, "[", "\\[");
76  boost::replace_all(pattern, "]", "\\]");
77  boost::replace_all(pattern, "*", "\\*");
78  boost::replace_all(pattern, "+", "\\+");
79  boost::replace_all(pattern, "?", "\\?");
80  boost::replace_all(pattern, "/", "\\/");
81 
82  boost::replace_all(pattern, "\\?", ".");
83  boost::replace_all(pattern, "\\*", ".*");
84 
85  boost::regex bpattern(pattern);
86 
87  return regex_match(text, bpattern);
88 }
89 
91 {
92  B2INFO("DQMHistAnalysisInputRootFile: beginRun called. Run: " << m_runList[m_run_idx]);
93  clearHistList();
94 }
95 
97 {
98  B2INFO("DQMHistAnalysisInputRootFile: event called.");
99 
100  sleep(m_interval);
101 
102  if (m_count > m_eventsList[m_run_idx]) {
103  m_run_idx++;
104  if (m_run_idx == m_runList.size()) {
105  m_eventMetaDataPtr.create();
106  m_eventMetaDataPtr->setEndOfData();
107  return;
108  }
109  m_count = 0;
110  if (m_file != nullptr) {
111  m_file->Close();
112  delete m_file;
113  }
114  m_file = new TFile(m_fileList[m_run_idx].c_str());
115  }
116 
117  // Clear only after EndOfRun check, otherwise we wont have any histograms for MiraBelle
118  // which expects analysis run in endRun function
120 
121  if (m_nullHistoMode) {
122  m_eventMetaDataPtr.create();
123  m_eventMetaDataPtr->setExperiment(m_expno);
125  m_eventMetaDataPtr->setEvent(m_count);
126  m_eventMetaDataPtr->setTime(0);
127  //setExpNr(m_expno); // redundant access from MetaData
128  //setRunNr(m_runno); // redundant access from MetaData
130  //ExtractRunType();
132  //ExtractEvent();
133 
134  B2INFO("DQMHistAnalysisInputRootFile: event finished. count: " << m_count);
135  m_count++;
136  return;
137  }
138 
139  std::vector<TH1*> hs;
140  unsigned long long int ts = 0;
141  m_file->cd();
142  TIter next(m_file->GetListOfKeys());
143  TKey* key = NULL;
144  while ((key = (TKey*)next())) {
145  TClass* cl = gROOT->GetClass(key->GetClassName());
146  if (ts == 0) ts = key->GetDatime().Convert();
147  if (!cl->InheritsFrom("TDirectory")) continue;
148  TDirectory* d = (TDirectory*)key->ReadObj();
149  std::string dirname = d->GetName();
150 
151  d->cd();
152  TIter nextd(d->GetListOfKeys());
153 
154  TKey* dkey;
155  while ((dkey = (TKey*)nextd())) {
156  TClass* dcl = gROOT->GetClass(dkey->GetClassName());
157  if (!dcl->InheritsFrom("TH1")) continue;
158  TH1* h = (TH1*)dkey->ReadObj();
159  if (h->InheritsFrom("TH2")) h->SetOption("col");
160  else h->SetOption("hist");
161  Double_t scale = 1.0 * m_count / m_eventsList[m_run_idx];
162  h->Scale(scale);
163  std::string hname = h->GetName();
164 
165  bool hpass = false;
166  if (m_histograms.size() == 0) {
167  hpass = true;
168  } else {
169  for (auto& hpattern : m_histograms) {
170  if (hnamePatternMatch(hpattern, dirname + "/" + hname)) {
171  hpass = true;
172  break;
173  }
174  }
175  }
176  if (!hpass) continue;
177 
178  if (hname.find("/") == std::string::npos) h->SetName((dirname + "/" + hname).c_str());
179  hs.push_back(h);
180  }
181  m_file->cd();
182  }
183 
184  // if no histograms are found in the sub-directories
185  // search the top folder
186  if (hs.size() == 0) {
187  TIter nexth(m_file->GetListOfKeys());
188  TKey* keyh = NULL;
189  while ((keyh = (TKey*)nexth())) {
190  TClass* cl = gROOT->GetClass(keyh->GetClassName());
191  TH1* h;
192  if (!cl->InheritsFrom("TH1")) continue;
193  h = (TH1*)keyh->ReadObj();
194 
195  bool hpass = false;
196  if (m_histograms.size() == 0) {
197  hpass = true;
198  } else {
199  for (auto& hpattern : m_histograms) {
200  if (hnamePatternMatch(hpattern, h->GetName())) {
201  hpass = true;
202  break;
203  }
204  }
205  }
206  if (!hpass) continue;
207 
208  hs.push_back(h);
209  Double_t scale = 1.0 * m_count / m_eventsList[m_run_idx];
210  h->Scale(scale);
211  }
212  }
213 
214  m_count++;
215  m_eventMetaDataPtr.create();
216  m_eventMetaDataPtr->setExperiment(m_expno);
218  m_eventMetaDataPtr->setEvent(m_count);
219  m_eventMetaDataPtr->setTime(ts * 1e9);
220 
221  //setExpNr(m_expno); // redundant access from MetaData
222  //setRunNr(m_runno); // redundant access from MetaData
223  if (m_runType == "") ExtractRunType(hs);
224  if (m_fillEvent <= 0) ExtractEvent(hs);
225 
226  // this code must be run after "event processed" has been extracted
227  for (size_t i = 0; i < hs.size(); i++) {
228  TH1* h = hs[i];
229  addHist("", h->GetName(), h);
230  B2DEBUG(1, "Found : " << h->GetName() << " : " << h->GetEntries());
231  }
232 
233  B2INFO("DQMHistAnalysisInputRootFile: event finished. count: " << m_count);
234 }
void initialize() override final
Initialize the module.
std::vector< unsigned int > m_eventsList
List of total number of events for each run.
std::vector< std::string > m_fileList
The list of names of the input root file.
bool m_nullHistoMode
Test mode for null histograms.
std::vector< std::string > m_histograms
List of histogram name patterns to process.
bool hnamePatternMatch(std::string pattern, std::string text)
Pattern match for histogram name.
TFile * m_file
The TFile object for the input file.
unsigned int m_interval
Time between two events in second.
void beginRun() override final
Called when entering a new run.
std::vector< unsigned int > m_runList
List of runs.
StoreObjPtr< EventMetaData > m_eventMetaDataPtr
Global EventMetaData for run number and event number.
unsigned int m_run_idx
Index in the list of runs, events and files.
The base class for the histogram analysis module.
static void clearHistList(void)
Clears the list of histograms.
void setRunType(std::string &t)
Set the Run Type.
void setEventProcessed(int e)
Set the number of processed events.
static bool addHist(const std::string &dirname, const std::string &histname, TH1 *h)
Add histogram.
void ExtractRunType(std::vector< TH1 * > &hs)
Extract Run Type from histogram title, called from input module.
static void initHistListBeforeEvent(void)
Reset the list of histograms.
void ExtractEvent(std::vector< TH1 * > &hs)
Extract event processed from daq histogram, called from input module.
REG_MODULE(arichBtest)
Register the Module.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
Abstract base class for different kinds of events.