Belle II Software  release-05-01-25
DQMHistAnalysisInput.cc
1 //+
2 // File : DQMHistAnalysisInput.cc
3 // Description :
4 //
5 // Author : Tomoyuki Konno, Tokyo Metropolitan Univerisity
6 // Date : 25 - Dec - 2015
7 //-
8 
9 
10 #include <dqm/analysis/modules/DQMHistAnalysisInput.h>
11 
12 #include <daq/slc/base/StringUtil.h>
13 #include <ctime>
14 
15 #include <TROOT.h>
16 #include <TKey.h>
17 
18 using namespace Belle2;
19 
20 //-----------------------------------------------------------------
21 // Register the Module
22 //-----------------------------------------------------------------
23 REG_MODULE(DQMHistAnalysisInput)
24 
25 //-----------------------------------------------------------------
26 // Implementation
27 //-----------------------------------------------------------------
28 
31 {
32  //Parameter definition
33  addParam("HistMemoryPath", m_mempath, "Path to Input Hist memory", std::string(""));
34  addParam("HistMemorySize", m_memsize, "Size of Input Hist memory", 10000000);
35  addParam("HistMemoryName", m_memname, "Name of Input Hist memory", std::string(""));
36  addParam("ShmId", m_shm_id, "ID of shared memory", 0);
37  addParam("SemId", m_sem_id, "ID of semaphore", 0);
38  addParam("RefreshInterval", m_interval, "Refresh interval of histograms", 10);
39  addParam("AutoCanvas", m_autocanvas, "Automatic creation of canvas", true);
40  addParam("AutoCanvasFolders", m_acfolders, "List of histograms to automatically create canvases, empty for all",
41  std::vector<std::string>());
42  addParam("ExcludeFolders", m_exclfolders, "List of folders to exclude from create canvases, empty for none, \"all\" for all",
43  std::vector<std::string>());
44  addParam("RemoveEmpty", m_remove_empty, "Remove empty histograms", false);
45  addParam("EnableRunInfo", m_enable_run_info, "Enable Run Info", false);
46  B2DEBUG(1, "DQMHistAnalysisInput: Constructor done.");
47 }
48 
49 
50 DQMHistAnalysisInputModule::~DQMHistAnalysisInputModule() { }
51 
53 {
54  if (m_memory != nullptr) delete m_memory;
55  if (m_mempath != "")
56  m_memory = new DqmMemFile(m_mempath.c_str());
57  else
60  m_c_info = new TCanvas("DQMInfo/c_info", "");
61  m_c_info->SetTitle("");
62  } else {
63  m_c_info = NULL;
64  }
65  m_eventMetaDataPtr.registerInDataStore();
66  B2INFO("DQMHistAnalysisInput: initialized.");
67 }
68 
69 
71 {
72  B2INFO("DQMHistAnalysisInput: beginRun called.");
73 }
74 
76 {
77  sleep(m_interval);
78  std::vector<TH1*> hs;
79  char mbstr[100];
80 
81  time_t now = time(0);
82  strftime(mbstr, sizeof(mbstr), "%c", localtime(&now));
83  B2INFO("[" << mbstr << "] before LoadMemFile");
84 
85  TMemFile* file = m_memory->LoadMemFile();
86 
87  now = time(0);
88  strftime(mbstr, sizeof(mbstr), "%c", localtime(&now));
89  B2INFO("[" << mbstr << "] after LoadMemFile");
90 
91  const TDatime& mt = file->GetModificationDate();
92  TDatime mmt(mt.Convert());
93  std::string expno("UNKNOWN"), runno("UNKNOWN"), rtype("UNKNOWN");
94 
95  file->cd();
96  TIter next(file->GetListOfKeys());
97  TKey* key = NULL;
98 
99  now = time(0);
100  strftime(mbstr, sizeof(mbstr), "%c", localtime(&now));
101  B2INFO("[" << mbstr << "] before input loop");
102 
103  while ((key = (TKey*)next())) {
104  TH1* h = (TH1*)key->ReadObj();
105  if (h == NULL) continue; // would be strange, but better check
106  if (m_remove_empty && h->GetEntries() == 0) continue;
107  // Remove ":" from folder name, workaround!
108  TString a = h->GetName();
109  a.ReplaceAll(":", "");
110  h->SetName(a);
111  B2DEBUG(1, "DQMHistAnalysisInput: get histo " << a.Data());
112 
113  if (StringUtil::split(a.Data(), '/').size() <= 1) continue;
114 
115  std::string dir_name = StringUtil::split(a.Data(), '/')[0];
116 
117  hs.push_back(h);
118  if (std::string(h->GetName()) == std::string("DQMInfo/expno")) expno = h->GetTitle();
119  if (std::string(h->GetName()) == std::string("DQMInfo/runno")) runno = h->GetTitle();
120  if (std::string(h->GetName()) == std::string("DQMInfo/rtype")) rtype = h->GetTitle();
121  if (m_autocanvas) {
122  StringList s = StringUtil::split(a.Data(), '/');
123 
124  bool give_canvas = false;
125  if (m_exclfolders.size() == 0) { //If none specified, canvases for all histograms
126  give_canvas = true;
127  } else {
128  bool in_excl_folder = false;
129  if (m_exclfolders.size() == 1 && m_exclfolders[0] == "all") {
130  in_excl_folder = true;
131  } else {
132  for (auto& excl_folder : m_exclfolders) {
133  if (excl_folder == s[0]) {
134  in_excl_folder = true;
135  break;
136  }
137  }
138  }
139 
140  if (in_excl_folder) {
141  for (auto& wanted_folder : m_acfolders) {
142  B2DEBUG(1, "==" << wanted_folder << "==" << s[0] << "==");
143  if (wanted_folder == std::string(h->GetName())) {
144  give_canvas = true;
145  break;
146  }
147  }
148  } else {
149  give_canvas = true;
150  }
151  }
152 
153  if (give_canvas) {
154  B2DEBUG(1, "Auto Hist->Canvas for " << a);
155  a.ReplaceAll("/", "_");
156  std::string name = a.Data();
157  if (m_cs.find(name) == m_cs.end()) {
158  if (s.size() > 1) {
159  std::string dirname = s[0];
160  std::string hname = s[1];
161  if ((dirname + "/" + hname) == "softwaretrigger/skim") hname = "skim_hlt";
162  TCanvas* c = new TCanvas((dirname + "/c_" + hname).c_str(), ("c_" + hname).c_str());
163  m_cs.insert(std::pair<std::string, TCanvas*>(name, c));
164  } else {
165  std::string hname = a.Data();
166  TCanvas* c = new TCanvas(("c_" + hname).c_str(), ("c_" + hname).c_str());
167  m_cs.insert(std::pair<std::string, TCanvas*>(name, c));
168  }
169  }
170  TCanvas* c = m_cs[name];
171  B2DEBUG(1, "DQMHistAnalysisInput: new canvas " << c->GetName());
172  c->cd();
173  if (h->GetDimension() == 1) {
174  if (h->GetMinimum() > 0) h->SetMinimum(0);
175  h->Draw("hist");
176  } else if (h->GetDimension() == 2) {
177  h->Draw("colz");
178  }
179  c->Update();
180  }
181  }
182  }
183 
184  now = time(0);
185  strftime(mbstr, sizeof(mbstr), "%c", localtime(&now));
186  B2INFO("[" << mbstr << "] after input loop");
187 
188  if (expno == std::string("UNKNOWN") || runno == std::string("UNKNOWN")) {
189  if (m_c_info != NULL) m_c_info->SetTitle((m_memname + ": Last Updated " + mmt.AsString()).c_str());
190  } else {
191  if (m_c_info != NULL) m_c_info->SetTitle((m_memname + ": Exp " + expno + ", Run " + runno + ", RunType " + rtype + ", Last Updated "
192  + mmt.AsString()).c_str());
193  m_expno = atoi(expno.c_str());
194  m_runno = atoi(runno.c_str());
195  }
196  B2INFO("DQMHistAnalysisInput: " << m_memname + ": Exp " + expno + ", Run " + runno + ", RunType " + rtype + ", Last Updated " +
197  mmt.AsString());
198 
199  resetHist();
200  for (size_t i = 0; i < hs.size(); i++) {
201  TH1* h = hs[i];
202  addHist("", h->GetName(), h);
203  B2DEBUG(1, "Found : " << h->GetName() << " : " << h->GetEntries());
204  std::string vname = h->GetName();
205  setFloatValue(vname + ".entries", h->GetEntries());
206  if (h->GetDimension() == 1) {
207  setFloatValue(vname + ".rms", h->GetRMS());
208  setFloatValue(vname + ".rmserr", h->GetRMSError());
209  setFloatValue(vname + ".mean", h->GetMean());
210  setFloatValue(vname + ".meanerr", h->GetMeanError());
211  } else if (h->GetDimension() == 2) {
212  setFloatValue(vname + ".xrms", h->GetRMS(1));
213  setFloatValue(vname + ".xrmserr", h->GetRMSError(1));
214  setFloatValue(vname + ".xmean", h->GetMean(1));
215  setFloatValue(vname + ".xmeanerr", h->GetMeanError(1));
216  setFloatValue(vname + ".yrms", h->GetRMS(2));
217  setFloatValue(vname + ".yrmserr", h->GetRMSError(2));
218  setFloatValue(vname + ".ymean", h->GetMean(2));
219  setFloatValue(vname + ".ymeanerr", h->GetMeanError(2));
220  }
221  }
222 
223  m_count++;
224  m_eventMetaDataPtr.create();
225  m_eventMetaDataPtr->setExperiment(m_expno);
226  m_eventMetaDataPtr->setRun(m_runno);
227  m_eventMetaDataPtr->setEvent(m_count);
228  m_eventMetaDataPtr->setTime(mt.Convert());
229 }
230 
232 {
233  B2INFO("DQMHistAnalysisInput : endRun called");
234 
235  TIter nextckey(gROOT->GetListOfCanvases());
236  TObject* cobj = NULL;
237 
238  while ((cobj = dynamic_cast<TObject*>(nextckey()))) {
239  if (cobj->IsA()->InheritsFrom("TCanvas")) {
240  (dynamic_cast<TCanvas*>(cobj))->Clear();
241  }
242  }
243 }
244 
245 
247 {
248  B2INFO("terminate called");
249 }
250 
Belle2::DQMHistAnalysisInputModule::m_remove_empty
bool m_remove_empty
Whether to remove empty histograms.
Definition: DQMHistAnalysisInput.h:68
Belle2::DQMHistAnalysisInputModule::m_memory
DqmMemFile * m_memory
Memory file to hold histograms.
Definition: DQMHistAnalysisInput.h:52
Belle2::DQMHistAnalysisInputModule::m_expno
unsigned int m_expno
Exp number.
Definition: DQMHistAnalysisInput.h:84
Belle2::DQMHistAnalysisModule::setFloatValue
static void setFloatValue(const std::string &parname, float vfloat)
Set the float value of the parameter.
Definition: DQMHistAnalysis.cc:181
Belle2::DQMHistAnalysisInputModule::m_exclfolders
std::vector< std::string > m_exclfolders
The list of folders which are excluded from automatically generate canvases.
Definition: DQMHistAnalysisInput.h:74
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::DQMHistAnalysisInputModule::m_sem_id
int m_sem_id
The semid for the shared memory.
Definition: DQMHistAnalysisInput.h:60
Belle2::DQMHistAnalysisInputModule::beginRun
virtual void beginRun() override
Module functions to be called from event process.
Definition: DQMHistAnalysisInput.cc:70
Belle2::DQMHistAnalysisInputModule::m_eventMetaDataPtr
StoreObjPtr< EventMetaData > m_eventMetaDataPtr
The metadata for each event.
Definition: DQMHistAnalysisInput.h:79
Belle2::DQMHistAnalysisInputModule::m_acfolders
std::vector< std::string > m_acfolders
The list of folders for which automatically generate canvases.
Definition: DQMHistAnalysisInput.h:72
Belle2::DQMHistAnalysisInputModule::m_c_info
TCanvas * m_c_info
The canvas hold the basic DQM info.
Definition: DQMHistAnalysisInput.h:76
Belle2::DQMHistAnalysisInputModule::m_memname
std::string m_memname
The name of the memory file (HLT or ExpressReco).
Definition: DQMHistAnalysisInput.h:56
Belle2::DQMHistAnalysisInputModule::initialize
virtual void initialize() override
Module functions to be called from main process.
Definition: DQMHistAnalysisInput.cc:52
Belle2::DQMHistAnalysisInputModule::m_mempath
std::string m_mempath
The name of the shared memory for the histograms.
Definition: DQMHistAnalysisInput.h:54
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisInputModule::m_interval
int m_interval
The refresh interval.
Definition: DQMHistAnalysisInput.h:64
Belle2::DQMHistAnalysisInputModule::m_shm_id
int m_shm_id
The shmid for the shared memory.
Definition: DQMHistAnalysisInput.h:58
Belle2::DQMHistAnalysisInputModule::m_count
unsigned int m_count
Event number.
Definition: DQMHistAnalysisInput.h:88
Belle2::DQMHistAnalysisInputModule::m_enable_run_info
bool m_enable_run_info
Whether to enable the run info to be displayed.
Definition: DQMHistAnalysisInput.h:70
Belle2::DQMHistAnalysisInputModule::m_runno
unsigned int m_runno
Run number.
Definition: DQMHistAnalysisInput.h:86
Belle2::DQMHistAnalysisInputModule::terminate
virtual void terminate() override
This method is called at the end of the event processing.
Definition: DQMHistAnalysisInput.cc:246
Belle2::DQMHistAnalysisModule::resetHist
static void resetHist()
Clear and reset the list of histograms.
Definition: DQMHistAnalysis.h:179
Belle2::DqmMemFile
Definition: DqmMemFile.h:28
Belle2::DQMHistAnalysisInputModule::m_autocanvas
bool m_autocanvas
Whether automatically generate canvases for histograms.
Definition: DQMHistAnalysisInput.h:66
Belle2::DQMHistAnalysisInputModule
Class definition for the output module of Sequential ROOT I/O.
Definition: DQMHistAnalysisInput.h:31
Belle2::DQMHistAnalysisInputModule::event
virtual void event() override
This method is the core of the module.
Definition: DQMHistAnalysisInput.cc:75
Belle2::DQMHistAnalysisModule::addHist
static void addHist(const std::string &dirname, const std::string &histname, TH1 *h)
Add histogram.
Definition: DQMHistAnalysis.cc:46
Belle2::DQMHistAnalysisInputModule::endRun
virtual void endRun() override
This method is called if the current run ends.
Definition: DQMHistAnalysisInput.cc:231
Belle2::DQMHistAnalysisInputModule::m_cs
std::map< std::string, TCanvas * > m_cs
The list of canvases for output.
Definition: DQMHistAnalysisInput.h:81
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27