Belle II Software development
DQMHistReferenceModule.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#include <dqm/analysis/modules/DQMHistReferenceModule.h>
10#include <TROOT.h>
11#include <TStyle.h>
12#include <TKey.h>
13
14using namespace std;
15using namespace Belle2;
16
17//-----------------------------------------------------------------
18// Register the Module
19//-----------------------------------------------------------------
20REG_MODULE(DQMHistReference);
21
22//-----------------------------------------------------------------
23// Implementation
24//-----------------------------------------------------------------
25
27{
28 //Parameter definition
29 addParam("ReferenceFile", m_referenceFileName, "Name of the reference histrogram files", string(""));
30 B2DEBUG(1, "DQMHistReference: Constructor done.");
31}
32
33
35
37{
38 B2DEBUG(1, "DQMHistReference: initialized.");
39}
40
42{
43 B2DEBUG(1, "DQMHistReference: beginRun called.");
44
46}
47
49{
50 TH1::AddDirectory(false); // do not store any histograms
51
52 B2DEBUG(1, "DQMHistReference: reading references from input root file");
53
54 string run_type = getRunType();
55 if (run_type == "") run_type = "default";
56
57 B2INFO("DQMHistReference: run_type " << run_type);
58
59 TFile* refFile = new TFile(m_referenceFileName.c_str(), "READ");
60
61 if (refFile->IsZombie()) {
62 B2INFO("DQMHistReference: reference file " << m_referenceFileName << " does not exist. No references will be used!");
63 refFile->Close();
64 delete refFile;
65 return;
66 }
67
68 B2INFO("DQMHistReference: use reference file " << m_referenceFileName);
69
70 TIter nextkey(refFile->GetListOfKeys());
71 TKey* key;
72 while ((key = (TKey*)nextkey())) {
73 if (key->IsFolder() && string(key->GetName()) == string("ref")) {
74 TDirectory* refdir = (TDirectory*)key->ReadObj(); // ReadObj -> I own it
75 TIter nextDetDir(refdir->GetListOfKeys());
76 TKey* detDir;
77 // detector folders
78 while ((detDir = (TKey*)nextDetDir())) {
79 if (!detDir->IsFolder()) continue;
80 TIter nextTypeDir(((TDirectory*)detDir->ReadObj())->GetListOfKeys());
81 TKey* typeDir;
82 TDirectory* foundDir = NULL;
83 // run type folders (get the run type corresponding folder or use default one)
84 while ((typeDir = (TKey*)nextTypeDir())) {
85 if (!typeDir->IsFolder()) continue;
86 if (string(typeDir->GetName()) == run_type) {
87 foundDir = (TDirectory*)typeDir->ReadObj(); // ReadObj -> I own it
88 break;
89 }
90 if (string(typeDir->GetName()) == "default") foundDir = (TDirectory*)typeDir->ReadObj(); // ReadObj -> I own it
91 }
92 string dirname = detDir->GetName();
93 if (!foundDir) {
94 B2INFO("No run type specific or default references available for " << dirname);
95 } else {
96 B2INFO("Reading reference histograms for " << dirname << " from run type folder: " << foundDir->GetName());
97
98 TIter next(foundDir->GetListOfKeys());
99 TKey* hh;
100
101 while ((hh = (TKey*)next())) {
102 if (hh->IsFolder()) continue;
103 if (gROOT->GetClass(key->GetClassName())->InheritsFrom("TH1")) {
104 addRefHist(dirname, (TH1*)key->ReadObj()); // ReadObj -> I own it, tranfer ownership to function;
105 }
106 }
107 delete foundDir; // always non-zero ... runtype or "default"
108 delete detDir; // always non-zero ... detector subdir name
109 }
110 }
111 delete refdir; // always non-zero ... "ref" folder
112 }
113 }
114
115 B2INFO("DQMHistReference: read references done");
116 refFile->Close();
117 delete refFile;
118}
119
121{
122 B2DEBUG(1, "DQMHistReference: event called");
123}
124
126{
127 B2DEBUG(1, "DQMHistReference: endRun called");
128}
129
131{
132 B2DEBUG(1, "DQMHistReference: terminate called");
133}
134
The base class for the histogram analysis module.
void addRefHist(const std::string &dirname, TH1 *hist)
Add reference histogram.
static const std::string & getRunType(void)
Get the list of the reference histograms.
void initialize() override final
Initializer.
void loadReferenceHistos()
Reads reference histograms from input root file.
void terminate() override final
This method is called at the end of the event processing.
void event() override final
This method is called for each event.
void endRun() override final
This method is called if the current run ends.
void beginRun() override final
Called when entering a new run.
std::string m_referenceFileName
Reference Histogram Root file name.
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:649
Abstract base class for different kinds of events.
STL namespace.