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