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