Belle II Software  release-05-02-19
DQMHistAnalysisOutputFile.cc
1 //+
2 // File : DQMHistAnalysisOutputFile.cc
3 // Description : DQM Analysis, dump histograms to file (as reference histograms)
4 //
5 // Author : B. Spruck
6 // Date : 25 - Mar - 2017
7 // based on work from Tomoyuki Konno, Tokyo Metropolitan Univerisity
8 //-
9 
10 
11 #include <dqm/analysis/modules/DQMHistAnalysisOutputFile.h>
12 #include <TROOT.h>
13 #include <TObject.h>
14 #include "TKey.h"
15 #include "TFile.h"
16 
17 using namespace std;
18 using namespace Belle2;
19 
20 //-----------------------------------------------------------------
21 // Register the Module
22 //-----------------------------------------------------------------
23 REG_MODULE(DQMHistAnalysisOutputFile)
24 
25 //-----------------------------------------------------------------
26 // Implementation
27 //-----------------------------------------------------------------
28 
31 {
32  //Parameter definition
33  addParam("HistoFile", m_filename, "Output Histogram Filename", std::string("histo.root"));
34  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
35  std::string("ref"));
36  addParam("SaveHistos", m_saveHistos, "Save Histos (default)", true);
37  addParam("SaveCanvases", m_saveCanvases, "Save Canvases (not default)", false);
38  B2DEBUG(20, "DQMHistAnalysisOutputFile: Constructor done.");
39 }
40 
41 
42 DQMHistAnalysisOutputFileModule::~DQMHistAnalysisOutputFileModule() { }
43 
44 void DQMHistAnalysisOutputFileModule::initialize()
45 {
46  B2DEBUG(20, "DQMHistAnalysisOutputFile: initialized.");
47 }
48 
49 
50 void DQMHistAnalysisOutputFileModule::beginRun()
51 {
52  B2DEBUG(20, "DQMHistAnalysisOutputFile: beginRun called.");
53 }
54 
55 
56 void DQMHistAnalysisOutputFileModule::event()
57 {
58  B2DEBUG(20, "DQMHistAnalysisOutputFile: event called.");
59 }
60 
61 void DQMHistAnalysisOutputFileModule::endRun()
62 {
63  B2INFO("DQMHistAnalysisOutputFile: endRun called");
64 
65  // Attention, we can not do that in Terminate, as then the memFile is already closed by previous task!
66  B2INFO("open file");
67  TFile f(m_filename.data(), "recreate");
68  if (f.IsOpen()) {
69  TDirectory* oldDir = gDirectory;
70  oldDir->mkdir(m_histogramDirectoryName.c_str());
71  oldDir->cd(m_histogramDirectoryName.c_str());
72 
73  if (m_saveCanvases) {
74  // we could loop over histos ... but i think the second one is better
75  TSeqCollection* seq;
76  seq = gROOT->GetListOfCanvases() ;
77  if (seq) {
78  B2INFO("found canvases");
79  TIter next(seq) ;
80  TObject* obj ;
81 
82  while ((obj = (TObject*)next())) {
83  if (obj->InheritsFrom("TCanvas")) {
84  B2INFO("Canvas name: " << obj->GetName() << " title " << obj->GetTitle());
85  obj->Write();
86  } else {
87  B2INFO("Others name: " << obj->GetName() << " title " << obj->GetTitle());
88  }
89  }
90  }
91  }
92  if (m_saveHistos) {
93  TSeqCollection* files;
94  files = gROOT->GetListOfFiles() ;
95  if (files) {
96  B2INFO("found keys");
97  TIter nextfile(files) ;
98  TObject* file ;
99 
100  while ((file = (TObject*)nextfile())) {
101  if (file->InheritsFrom("TFile")) {
102  B2INFO("File name: " << file->GetName() << " title " << file->GetTitle());
103  if (file == &f || file->GetName() == m_filename) continue;
104 
105  TList* list = ((TFile*)file)->GetListOfKeys() ;
106  if (list) {
107  TIter next(list) ;
108  TKey* key ;
109  TObject* obj ;
110 
111  while ((key = (TKey*)next())) {
112  TString skey(key->GetClassName());
113  if (skey.BeginsWith(TString("Belle2::"))) continue;
114  TClass clkey(key->GetClassName());
115  if (clkey.InheritsFrom("TH1")) {
116  obj = key->ReadObj() ;
117  B2INFO("Histo name: " << obj->GetName() << " title " << obj->GetTitle());
118  TDirectory* old, *d;
119  d = old = gDirectory;
120  TString myl = obj->GetName();
121  TString tok;
122 // workaround for changes in histmanager/Input module ... TODO remove if not needed anymore
123 // std::string a=myl.Data();
124 // a.erase(std::remove(a.begin(), a.end(), ':'), a.end());
125 // myl=a.data();
126  Ssiz_t from = 0;
127  while (myl.Tokenize(tok, from, "/")) {
128  TString dummy;
129  Ssiz_t fr;
130  fr = from;
131  if (myl.Tokenize(dummy, fr, "/")) { // check if its the last one
132  auto e = d->GetDirectory(tok);
133  if (e) {
134  d = e;
135  d->cd();
136  } else {
137  d->mkdir(tok);
138  d->cd(tok);
139  d = gDirectory;
140  }
141  } else {
142  break;
143  }
144  }
145  ((TH1*)obj)->SetName(tok);
146  obj->Write();
147  old->cd();
148  }
149  }
150  }
151  } else {
152  B2INFO("Others name: " << file->GetName() << " title " << file->GetTitle());
153  }
154  }
155  }
156  }
157 
158 
159  f.Write();
160  f.Close();
161  }
162 }
163 
164 
165 void DQMHistAnalysisOutputFileModule::terminate()
166 {
167  B2INFO("DQMHistAnalysisOutputFile: terminate called");
168 // Attention, we can not do that in Terminate, as then the memFile is already closed by previous task!
169 }
170 
Belle2::DQMHistAnalysisOutputFileModule
Class definition for the output module of Sequential ROOT I/O.
Definition: DQMHistAnalysisOutputFile.h:22
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
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