Belle II Software  release-05-01-25
DQMHistAnalysisOutputMonObj.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Luka Santelj *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 
12 #include <dqm/analysis/modules/DQMHistAnalysisOutputMonObj.h>
13 #include <framework/dataobjects/EventMetaData.h>
14 #include <framework/datastore/DataStore.h>
15 #include <framework/datastore/StoreObjPtr.h>
16 #include "TFile.h"
17 #include "TTree.h"
18 #include "TString.h"
19 #include <time.h>
20 
21 
22 using namespace std;
23 using namespace Belle2;
24 
25 //-----------------------------------------------------------------
26 // Register the Module
27 //-----------------------------------------------------------------
28 REG_MODULE(DQMHistAnalysisOutputMonObj)
29 
30 //-----------------------------------------------------------------
31 // Implementation
32 //-----------------------------------------------------------------
33 
36 {
37  //Parameter definition
38  addParam("Filename", m_filename, "Output root filename (if not set mon_e{exp}r{run}.root is used", std::string(""));
39  addParam("TreeFile", m_treeFile, "If set, entry to run summary TTree from TreeFile is made", std::string(""));
40  addParam("ProcID", m_procID, "Processing id (online,proc10, etc.)", std::string("online"));
41  addParam("run", m_run, "Run number", 0);
42  addParam("exp", m_exp, "Experiment number", 0);
43  addParam("nevt", m_nevt, "Number of events", 0);
44  addParam("runtype", m_runtype, "Run type", std::string(""));
45 
46  B2DEBUG(20, "DQMHistAnalysisOutputMonObj: Constructor done.");
47 }
48 
49 
50 DQMHistAnalysisOutputMonObjModule::~DQMHistAnalysisOutputMonObjModule() { }
51 
52 void DQMHistAnalysisOutputMonObjModule::initialize()
53 {
54  B2DEBUG(20, "DQMHistAnalysisOutputMonObj: initialized.");
55  // create file metadata
56  m_metaData = new DQMFileMetaData();
57  m_metaData->setProcessingID(m_procID);
58 }
59 
60 
61 void DQMHistAnalysisOutputMonObjModule::beginRun()
62 {
63  B2DEBUG(20, "DQMHistAnalysisOutputMonObj: beginRun called.");
64 }
65 
66 
67 void DQMHistAnalysisOutputMonObjModule::event()
68 {
69  B2DEBUG(20, "DQMHistAnalysisOutputMonObj: event called.");
70 }
71 
72 void DQMHistAnalysisOutputMonObjModule::endRun()
73 {
74  B2INFO("DQMHistAnalysisOutputMonObj: endRun called");
75 
76  StoreObjPtr<EventMetaData> lastEvtMeta;
77 
78  B2INFO("open file");
79  TH1* hrun = findHist("DQMInfo/runno");
80  TH1* hexp = findHist("DQMInfo/expno");
81 
82  int run = hrun ? std::stoi(hrun->GetTitle()) : m_run;
83  int exp = hexp ? std::stoi(hexp->GetTitle()) : m_exp;
84  TString fname;
85  if (m_filename.length()) fname = m_filename;
86  else fname = TString::Format("mon_e%04dr%06d_%s.root", exp, run, m_procID.c_str());
87 
88  TH1* runtype = findHist("DQMInfo/rtype");
89  if (runtype) m_metaData->setRunType(std::string(runtype->GetTitle()));
90  else m_metaData->setRunType(m_runtype);
91  TH1* hnevt = findHist("DAQ/Nevent");
92  if (hnevt) m_metaData->setNEvents(hnevt->GetEntries());
93  else m_metaData->setNEvents(m_nevt);
94 
95  TFile f(fname, "NEW");
96 
97  if (f.IsZombie()) {
98  B2WARNING("File " << LogVar("MonitoringObject file",
99  fname) << " already exists and it will not be rewritten. If desired please delete file and re-run.");
100  return;
101  }
102 
103  // set meta data info
104  // m_metaData->setNEvents(lastEvtMeta->getEvent());
105  m_metaData->setExperimentRun(exp, run);
106  time_t ts = lastEvtMeta->getTime() / 1e9;
107  struct tm* timeinfo;
108  timeinfo = localtime(&ts);
109  m_metaData->setRunDate(asctime(timeinfo));
110 
111  m_metaData->Write();
112  // get list of existing monitoring objects
113  const MonObjList& objts = getMonObjList();
114  // write them to the output file
115  for (const auto& obj : objts)(obj.second)->Write();
116 
117  f.Close();
118 
119  if (m_treeFile.length() > 0) addTreeEntry();
120 
121 }
122 
123 void DQMHistAnalysisOutputMonObjModule::addTreeEntry()
124 {
125 
126  TFile* treeFile = new TFile(m_treeFile.c_str(), "update");
127  auto* tree = (TTree*)treeFile->Get("tree");
128 
129  if (tree == NULL) tree = new TTree("tree", "tree");
130 
131  int run = m_metaData->getRun();
132  int expe = m_metaData->getExperiment();
133  int nevt = m_metaData->getNEvents();
134  //int rune = 0;
135  //int expee = 0;
136  char* rel = const_cast<char*>(m_metaData->getRelease().c_str());
137  char* db = const_cast<char*>(m_metaData->getDatabaseGlobalTag().c_str());
138  char* datee = const_cast<char*>(m_metaData->getRunDate().c_str());
139  char* rtype = const_cast<char*>(m_metaData->getRunType().c_str());
140  char* procID = const_cast<char*>(m_metaData->getProcessingID().c_str());
141 
142  auto b_run = tree->GetBranch("run");
143  auto b_exp = tree->GetBranch("exp");
144  auto b_release = tree->GetBranch("release");
145  auto b_gt = tree->GetBranch("gt");
146  auto b_datetime = tree->GetBranch("datetime");
147  auto b_rtype = tree->GetBranch("rtype");
148  auto b_procID = tree->GetBranch("procID");
149  auto b_nevt = tree->GetBranch("nevt");
150 
151  // this still needs to be sorted out
152  /*if(b_run){
153  b_run->SetAddress(&rune);
154  b_exp->SetAddress(&expee);
155  for(int ie = 0; ie<b_run->GetEntries(); ie++){
156  b_run->GetEntry(ir);
157  b_exp->GetEntry(ir);
158  if(rune == run && expee = expe){
159 
160  }
161  }
162  }*/
163 
164 
165  if (!b_run) tree->Branch("run", &run, "run/I");
166  else b_run->SetAddress(&run);
167  if (!b_exp) tree->Branch("exp", &expe, "exp/I");
168  else b_exp->SetAddress(&expe);
169  if (!b_nevt) tree->Branch("nevt", &nevt, "nevt/I");
170  else b_nevt->SetAddress(&nevt);
171  if (!b_release) tree->Branch("release", rel, "release/C");
172  else b_release->SetAddress(rel);
173  if (!b_gt) tree->Branch("gt", db, "gt/C");
174  else b_gt->SetAddress(db);
175  if (!b_datetime) tree->Branch("datetime", datee, "datetime/C");
176  else b_datetime->SetAddress(datee);
177  if (!b_rtype) tree->Branch("rtype", rtype, "rtype/C");
178  else b_rtype->SetAddress(rtype);
179  if (!b_procID) tree->Branch("procID", procID, "procID/C");
180  else b_procID->SetAddress(procID);
181 
182 
183  const MonObjList& objts = getMonObjList();
184  // write them to the output file
185  for (const auto& obj : objts) {
186  std::map<std::string, float>& vars = const_cast<std::map<std::string, float>&>((obj.second)->getVariables());
187  std::map<std::string, float>& upErr = const_cast<std::map<std::string, float>&>((obj.second)->getUpError());
188  std::map<std::string, float>& lowErr = const_cast<std::map<std::string, float>&>((obj.second)->getLowError());
189 
190  const std::vector<std::pair<std::string, std::string>>& strVars = (obj.second)->getStringVariables();
191 
192  for (auto& var : vars) {
193  std::string brname = obj.first + "_" + var.first;
194  auto branch = tree->GetBranch((brname).c_str());
195  if (!branch) {
196  branch = tree->Branch((brname).c_str(), &(var.second));
197  fillBranch(branch);
198  } else branch->SetAddress(&(var.second));
199 
200  auto vvE1 = upErr.find(var.first);
201  auto vvE2 = lowErr.find(var.first);
202 
203  if (vvE1 != upErr.end() && vvE2 == lowErr.end()) {
204  auto errBranch = tree->GetBranch((brname).c_str() + TString("_err"));
205  if (!errBranch) {
206  errBranch = tree->Branch((brname).c_str() + TString("_err"), &(vvE1->second));
207  fillBranch(errBranch);
208  } else errBranch->SetAddress(&(vvE1->second));
209  }
210 
211  if (vvE1 != upErr.end() && vvE2 != lowErr.end()) {
212  auto errBranch1 = tree->GetBranch((brname).c_str() + TString("_upErr"));
213  if (!errBranch1) {
214  errBranch1 = tree->Branch((brname).c_str() + TString("_upErr"), &(vvE1->second));
215  fillBranch(errBranch1);
216  } else errBranch1->SetAddress(&(vvE1->second));
217 
218  auto errBranch2 = tree->GetBranch((brname).c_str() + TString("_dwErr"));
219  if (!errBranch2) {
220  errBranch2 = tree->Branch((brname).c_str() + TString("_dwErr"), &(vvE2->second));
221  fillBranch(errBranch2);
222  } else errBranch2->SetAddress(&(vvE2->second));
223 
224  }
225  }
226 
227  for (auto& var : strVars) {
228  std::string brname = obj.first + "_" + var.first;
229  char* cc = const_cast<char*>((var.second).c_str());
230  auto branch = tree->GetBranch((brname).c_str());
231  if (!branch) {
232  std::string ty = brname + "/C";
233  branch = tree->Branch((brname).c_str(), cc, ty.c_str());
234  fillBranch(branch);
235  } else branch->SetAddress(cc);
236  }
237  }
238 
239  tree->Fill();
240  tree->Write(0, TObject::kWriteDelete, 0);
241  treeFile->Close();
242 
243 }
244 
245 void DQMHistAnalysisOutputMonObjModule::fillBranch(TBranch* branch)
246 {
247  TTree* tree = (TTree*)branch->GetTree();
248  int nentr = tree->GetEntries();
249  for (int i = 0; i < nentr; i++) {
250  tree->GetEntry(i);
251  branch->Fill();
252  }
253 }
254 
255 
256 
257 void DQMHistAnalysisOutputMonObjModule::terminate()
258 {
259  B2INFO("DQMHistAnalysisOutputMonObj: terminate called");
260 // Attention, we can not do that in Terminate, as then the memFile is already closed by previous task!
261 
262 }
263 
Belle2::DQMHistAnalysisOutputMonObjModule
Class definition for the module to store MonitoringObject to output root file.
Definition: DQMHistAnalysisOutputMonObj.h:33
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::DQMFileMetaData
Metadata information about a DQM file.
Definition: DQMFileMetaData.h:37
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
LogVar
Class to store variables with their name which were sent to the logging service.
Definition: LogVariableStream.h:24
Belle2::DQMHistAnalysisModule::MonObjList
std::map< std::string, MonitoringObject * > MonObjList
The type of list of MonitoringObjects.
Definition: DQMHistAnalysis.h:70
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27