Belle II Software  release-08-01-10
DQMHistAnalysisPXDInjection.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 : DQMHistAnalysisPXDInjection.cc
10 // Description : DQM module, which gives histograms showing occupancies after injection
11 //-
12 
13 
14 #include <dqm/analysis/modules/DQMHistAnalysisPXDInjection.h>
15 #include <TROOT.h>
16 
17 using namespace std;
18 using namespace Belle2;
19 
20 //-----------------------------------------------------------------
21 // Register the Module
22 //-----------------------------------------------------------------
23 REG_MODULE(DQMHistAnalysisPXDInjection);
24 
25 //-----------------------------------------------------------------
26 // Implementation
27 //-----------------------------------------------------------------
28 
29 DQMHistAnalysisPXDInjectionModule::DQMHistAnalysisPXDInjectionModule() : DQMHistAnalysisModule()
30 {
31  // This module CAN NOT be run in parallel!
32  setDescription("DQM Analysis for PXD Hits after Injection");
33 
34  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms were placed",
35  std::string("PXDINJ"));
36  B2DEBUG(1, "DQMHistAnalysisPXDInjection: Constructor done.");
37 }
38 
40 {
41 
42  gROOT->cd(); // this seems to be important, or strange things happen
43 
45  TString("/c_PXDInjectionLER")); // different name, as we have the other module
47  TString("/c_PXDInjectionHER")); // which makes this hist for ALL detectors (SVD, ECl, TOP)
48  m_hInjectionLERPXD = new TH1F("HitInjectionLERPXD", "PXD Hits after LER Injection;Time in #mus;Mean Hits/event", 4000, 0, 20000);
49  m_hInjectionHERPXD = new TH1F("HitInjectionHERPXD", "PXD Hits after HER Injection;Time in #mus;Mean Hits/event", 4000, 0, 20000);
50  const VXD::GeoCache& vxdGeometry = VXD::GeoCache::getInstance();
51  std::vector<VxdID> vxdsensors = vxdGeometry.getListOfSensors();
52  for (VxdID& avxdid : vxdsensors) {
53  VXD::SensorInfoBase info = vxdGeometry.getSensorInfo(avxdid);
54  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
55  m_sensors.push_back(avxdid);
56  }
57 
58  for (VxdID& avxdid : m_sensors) {
59  TString buff = (std::string)avxdid;
60  TString bufful = buff;
61  bufful.ReplaceAll(".", "_");
62 
63  m_cInjectionLERPXDMod[avxdid] = new TCanvas(m_histogramDirectoryName + "/c_PXDInjectionLER_" + bufful);
64  m_cInjectionHERPXDMod[avxdid] = new TCanvas(m_histogramDirectoryName + "/c_PXDInjectionHER_" + bufful);
65  m_hInjectionLERPXDMod[avxdid] = new TH1F("HitInjectionLERPXD_" + bufful,
66  "PXD Hits after LER Injection " + buff + "/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
67  m_hInjectionHERPXDMod[avxdid] = new TH1F("HitInjectionHERPXD_" + bufful,
68  "PXD Hits after HER Injection " + buff + "/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
69  m_cInjectionLERPXDModNorm[avxdid] = new TCanvas(m_histogramDirectoryName + "/c_PXDInjectionLERNorm_" + bufful);
70  m_cInjectionHERPXDModNorm[avxdid] = new TCanvas(m_histogramDirectoryName + "/c_PXDInjectionHERNorm_" + bufful);
71  m_hInjectionLERPXDModNorm[avxdid] = new TH1F("HitInjectionLERPXDNorm_" + bufful,
72  "PXD Hits after LER Injection " + buff + " normalized to 1.1.x/Time;Time in #mus;factor", 4000, 0, 20000);
73  m_hInjectionHERPXDModNorm[avxdid] = new TH1F("HitInjectionHERPXDNorm_" + bufful,
74  "PXD Hits after HER Injection " + buff + " normalized to 1.1.x/Time;Time in #mus;factor", 4000, 0, 20000);
75  }
76 
77  B2DEBUG(1, "DQMHistAnalysisPXDInjection: initialized.");
78 }
79 
80 
82 {
83  B2DEBUG(1, "DQMHistAnalysisPXDInjection: beginRun called.");
84 
85 // m_cInjectionLERPXD->Clear(); // FIXME, unclear if this lets to crashes on new run?
86 // m_cInjectionLERPXDOcc->Clear();
87 // m_cInjectionLERECL->Clear();
88 // m_cInjectionHERPXD->Clear();
89 // m_cInjectionHERPXDOcc->Clear();
90 // m_cInjectionHERECL->Clear();
91 }
92 
93 
95 {
96  TH1* Triggers = nullptr;
97  // cppcheck-suppress unreadVariable
98  TString locationHits = "";
99  TString locationTriggers = "PXDEOccInjLER";
100  if (m_histogramDirectoryName != "") {
101  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
102  }
103  Triggers = (TH1*)findHist(locationTriggers.Data());
104 
105  //Finding only one of them should only happen in very strange situations...
106  //m_nodes[0].histo = Triggers;
107  if (Triggers) {
108  TH1* Hits = nullptr, *RefMod_fw = nullptr, *RefMod_bw = nullptr;
109  locationHits = "PXDOccInjLER";
110  if (m_histogramDirectoryName != "") {
111  locationHits = m_histogramDirectoryName + "/" + locationHits;
112  }
113  Hits = (TH1*)findHist(locationHits.Data());
114  if (Hits) {
115  m_hInjectionLERPXD->Divide(Hits, Triggers);
116  }
117  locationHits = "PXDOccInjLER_1_1_1";
118  if (m_histogramDirectoryName != "") {
119  locationHits = m_histogramDirectoryName + "/" + locationHits;
120  }
121  RefMod_fw = (TH1*)findHist(locationHits.Data());
122  locationHits = "PXDOccInjLER_1_1_2";
123  if (m_histogramDirectoryName != "") {
124  locationHits = m_histogramDirectoryName + "/" + locationHits;
125  }
126  RefMod_bw = (TH1*)findHist(locationHits.Data());
127 
128  // assume trigger is the same for all modules (makes sense :-)
129  for (VxdID& avxdid : m_sensors) {
130  TString buff = (std::string)avxdid;
131  TString bufful = buff;
132  bufful.ReplaceAll(".", "_");
133 
134  locationHits = "PXDOccInjLER_" + bufful;
135  if (m_histogramDirectoryName != "") {
136  locationHits = m_histogramDirectoryName + "/" + locationHits;
137  }
138  Hits = (TH1*)findHist(locationHits.Data());
139  if (Hits) {
140  m_hInjectionLERPXDMod[avxdid]->Divide(Hits, Triggers);
141  if (avxdid.getSensorNumber() == 1 && RefMod_fw) {
142  m_hInjectionLERPXDModNorm[avxdid]->Divide(Hits, RefMod_fw);
143  } else if (avxdid.getSensorNumber() == 2 && RefMod_bw) {
144  m_hInjectionLERPXDModNorm[avxdid]->Divide(Hits, RefMod_bw);
145  }
146  }
147 
148  }
149 
150  }
151 
152  locationTriggers = "PXDEOccInjHER";
153  if (m_histogramDirectoryName != "") {
154  locationTriggers = m_histogramDirectoryName + "/" + locationTriggers;
155  }
156  Triggers = (TH1*)findHist(locationTriggers.Data());
157 
158  //Finding only one of them should only happen in very strange situations...
159  //m_nodes[3].histo = Triggers;
160  if (Triggers) {
161  TH1* Hits = nullptr, *RefMod_fw = nullptr, *RefMod_bw = nullptr;
162  locationHits = "PXDOccInjHER";
163  if (m_histogramDirectoryName != "") {
164  locationHits = m_histogramDirectoryName + "/" + locationHits;
165  }
166  Hits = (TH1*)findHist(locationHits.Data());
167  if (Hits) {
168  m_hInjectionHERPXD->Divide(Hits, Triggers);
169  }
170  locationHits = "PXDOccInjHER_1_1_1";
171  if (m_histogramDirectoryName != "") {
172  locationHits = m_histogramDirectoryName + "/" + locationHits;
173  }
174  RefMod_fw = (TH1*)findHist(locationHits.Data());
175  locationHits = "PXDOccInjHER_1_1_2";
176  if (m_histogramDirectoryName != "") {
177  locationHits = m_histogramDirectoryName + "/" + locationHits;
178  }
179  RefMod_bw = (TH1*)findHist(locationHits.Data());
180  // assume trigger is the same for all modules (makes sense :-)
181  for (VxdID& avxdid : m_sensors) {
182  TString buff = (std::string)avxdid;
183  TString bufful = buff;
184  bufful.ReplaceAll(".", "_");
185  locationHits = "PXDOccInjHER_" + bufful;
186  if (m_histogramDirectoryName != "") {
187  locationHits = m_histogramDirectoryName + "/" + locationHits;
188  }
189  Hits = (TH1*)findHist(locationHits.Data());
190  if (Hits) {
191  m_hInjectionHERPXDMod[avxdid]->Divide(Hits, Triggers);
192  if (avxdid.getSensorNumber() == 1 && RefMod_fw) {
193  m_hInjectionHERPXDModNorm[avxdid]->Divide(Hits, RefMod_fw);
194  } else if (avxdid.getSensorNumber() == 2 && RefMod_bw) {
195  m_hInjectionHERPXDModNorm[avxdid]->Divide(Hits, RefMod_bw);
196  }
197  }
198 
199 
200  }
201  }
202 
203  m_cInjectionLERPXD->Clear();
204  m_cInjectionLERPXD->cd(0);
205  m_cInjectionLERPXD->Pad()->SetLogy();
206  m_hInjectionLERPXD->Draw("hist");
207 
208  m_cInjectionHERPXD->Clear();
209  m_cInjectionHERPXD->cd(0);
210  m_cInjectionHERPXD->Pad()->SetLogy();
211  m_hInjectionHERPXD->Draw("hist");
212 
213  for (VxdID& avxdid : m_sensors) {
214  m_cInjectionHERPXDMod[avxdid]->Clear();
215  m_cInjectionHERPXDMod[avxdid]->cd(0);
216  m_cInjectionHERPXDMod[avxdid]->Pad()->SetLogy();
217  m_hInjectionHERPXDMod[avxdid]->Draw("hist");
218  m_cInjectionHERPXDModNorm[avxdid]->Clear();
219  m_cInjectionHERPXDModNorm[avxdid]->cd(0);
220  m_hInjectionHERPXDModNorm[avxdid]->Draw("hist");
221  m_cInjectionLERPXDMod[avxdid]->Clear();
222  m_cInjectionLERPXDMod[avxdid]->cd(0);
223  m_cInjectionLERPXDMod[avxdid]->Pad()->SetLogy();
224  m_hInjectionLERPXDMod[avxdid]->Draw("hist");
225  m_cInjectionLERPXDModNorm[avxdid]->Clear();
226  m_cInjectionLERPXDModNorm[avxdid]->cd(0);
227  m_hInjectionLERPXDModNorm[avxdid]->Draw("hist");
228  }
229 }
230 
The base class for the histogram analysis module.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
std::map< VxdID, TCanvas * > m_cInjectionLERPXDModNorm
Canvases per sensor for LER normalized.
std::map< VxdID, TH1F * > m_hInjectionHERPXDModNorm
Histogram per sensor for HER normalized.
void initialize(void) override final
Initializer.
std::map< VxdID, TCanvas * > m_cInjectionHERPXDModNorm
Canvases per sensor for HER normalized.
std::string m_histogramDirectoryName
name of histogram directory
std::map< VxdID, TCanvas * > m_cInjectionLERPXDMod
Canvases per sensor for LER.
std::map< VxdID, TH1F * > m_hInjectionHERPXDMod
Histogram per sensor for HER.
std::map< VxdID, TH1F * > m_hInjectionLERPXDMod
Histogram per sensor for LER.
std::map< VxdID, TH1F * > m_hInjectionLERPXDModNorm
Histogram per sensor for LER normalized.
std::vector< VxdID > m_sensors
List of PXD sensors.
std::map< VxdID, TCanvas * > m_cInjectionHERPXDMod
Canvases per sensor for HER.
void beginRun(void) override final
Called when entering a new run.
void event(void) override final
This method is called for each event.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition: GeoCache.cc:59
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
Base class to provide Sensor Information for PXD and SVD.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
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
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.