Belle II Software  release-08-01-10
PXDInjectionDQMModule.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 <pxd/modules/pxdDQM/PXDInjectionDQMModule.h>
10 #include "TDirectory.h"
11 
12 using namespace std;
13 using namespace Belle2;
14 using namespace Belle2::PXD;
15 using namespace Belle2::VXD;
16 
17 //-----------------------------------------------------------------
18 // Register the Module
19 //-----------------------------------------------------------------
20 REG_MODULE(PXDInjectionDQM);
21 
22 //-----------------------------------------------------------------
23 // Implementation
24 //-----------------------------------------------------------------
25 
26 PXDInjectionDQMModule::PXDInjectionDQMModule() : HistoModule(), m_vxdGeometry(VXD::GeoCache::getInstance())
27 {
28  //Set module properties
29  setDescription("Monitor Occupancy after Injection");
31  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
32  std::string("PXDINJ"));
33  addParam("PXDRawHitsName", m_PXDRawHitsName, "Name of PXD raw hits", std::string(""));
34  addParam("PXDClusterName", m_PXDClustersName, "Name of PXD clusters", std::string(""));
35  addParam("eachModule", m_eachModule, "create for each module", false);
36  addParam("offlineStudy", m_offlineStudy, "use finest binning and larger range", false);
37  addParam("useClusters", m_useClusters, "use cluster instead of raw hits", false);
38  addParam("createMaxHist", m_createMaxHist, "create histo with max occupancy (not mp save!!!)", false);
39  addParam("createGateHist", m_createGateHist, "create 2d histo with gate against occupancy", false);
40 
41 }
42 
44 {
45  TDirectory* oldDir = gDirectory;
46  if (m_histogramDirectoryName != "") {
47  oldDir->mkdir(m_histogramDirectoryName.c_str());// do not rely on return value, might be ZERO
48  oldDir->cd(m_histogramDirectoryName.c_str());//changing to the right directory
49  }
50 
51  if (m_offlineStudy) {
52  hOccAfterInjLER = new TH1F("PXDOccInjLER", "PXDOccInjLER/Time;Time in #mus;Count/Time (0.5 #mus bins)", 100000, 0, 50000);
53  hOccAfterInjHER = new TH1F("PXDOccInjHER", "PXDOccInjHER/Time;Time in #mus;Count/Time (0.5 #mus bins)", 100000, 0, 50000);
54  hEOccAfterInjLER = new TH1I("PXDEOccInjLER", "PXDEOccInjLER/Time;Time in #mus;Triggers/Time (0.5 #mus bins)", 100000, 0, 50000);
55  hEOccAfterInjHER = new TH1I("PXDEOccInjHER", "PXDEOccInjHER/Time;Time in #mus;Triggers/Time (0.5 #mus bins)", 100000, 0, 50000);
56  if (m_createMaxHist) {
57  hMaxOccAfterInjLER = new TH1F("PXDMaxOccInjLER", "PXDMaxOccInjLER/Time;Time in #mus;Triggers/Time (0.5 #mus bins)", 100000, 0,
58  50000);
59  hMaxOccAfterInjHER = new TH1F("PXDMaxOccInjHER", "PXDMaxOccInjHER/Time;Time in #mus;Triggers/Time (0.5 #mus bins)", 100000, 0,
60  50000);
61  }
62  if (m_createGateHist) {
63  hOccAfterInjLERGate = new TH2F("PXDOccInjLERGate", "PXDOccInjLERGate;Time;Gate", 1000, 0, 10000, 192, 0, 192);
64  hOccAfterInjHERGate = new TH2F("PXDOccInjHERGate", "PXDOccInjHERGate;Time;Gate", 1000, 0, 10000, 192, 0, 192);
65  }
66  } else {
67  hOccAfterInjLER = new TH1F("PXDOccInjLER", "PXDOccInjLER/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
68  hOccAfterInjHER = new TH1F("PXDOccInjHER", "PXDOccInjHER/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
69  hEOccAfterInjLER = new TH1I("PXDEOccInjLER", "PXDEOccInjLER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0, 20000);
70  hEOccAfterInjHER = new TH1I("PXDEOccInjHER", "PXDEOccInjHER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0, 20000);
71  if (m_createMaxHist) {
72  hMaxOccAfterInjLER = new TH1F("PXDMaxOccInjLER", "PXDMaxOccInjLER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0, 20000);
73  hMaxOccAfterInjHER = new TH1F("PXDMaxOccInjHER", "PXDMaxOccInjHER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0, 20000);
74  }
75  if (m_createGateHist) {
76  hOccAfterInjLERGate = new TH2F("PXDOccInjLERGate", "PXDOccInjLERGate;Time;Gate", 1000, 0, 10000, 192, 0, 192);
77  hOccAfterInjHERGate = new TH2F("PXDOccInjHERGate", "PXDOccInjHERGate;Time;Gate", 1000, 0, 10000, 192, 0, 192);
78  }
79  }
80 
81  if (m_eachModule) {
82  std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
83  for (VxdID& avxdid : sensors) {
85  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
86  // Only interested in PXD sensors
87 
88  TString buff = (std::string)avxdid;
89  TString bufful = buff;
90  bufful.ReplaceAll(".", "_");
91 
94  if (m_offlineStudy) {
95  hOccModAfterInjLER[avxdid] = new TH1F("PXDOccInjLER_" + bufful,
96  "PXDOccModInjLER " + buff + "/Time;Time in #mus;Count/Time (0.5 #mus bins)", 100000, 0, 50000);
97  hOccModAfterInjHER[avxdid] = new TH1F("PXDOccInjHER_" + bufful,
98  "PXDOccModInjHER " + buff + "/Time;Time in #mus;Count/Time (0.5 #mus bins)", 100000, 0, 50000);
99  if (m_createMaxHist) {
100  hMaxOccModAfterInjLER[avxdid] = new TH1F("PXDMaxOccInjLER_" + bufful,
101  "PXDMaxOccModInjLER " + buff + "/Time;Time in #mus;Count/Time (0.5 #mus bins)", 100000, 0, 50000);
102  hMaxOccModAfterInjHER[avxdid] = new TH1F("PXDMaxOccInjHER_" + bufful,
103  "PXDMaxOccModInjHER " + buff + "/Time;Time in #mus;Count/Time (0.5 #mus bins)", 100000, 0, 50000);
104  }
105  if (m_createGateHist) {
106  hOccModAfterInjLERGate[avxdid] = new TH2F("PXDOccInjLERGate_" + bufful, "PXDOccInjLERGate " + buff + ";Time;Gate", 1000, 0, 10000,
107  192, 0,
108  192);
109  hOccModAfterInjHERGate[avxdid] = new TH2F("PXDOccInjHERGate_" + bufful, "PXDOccInjHERGate " + buff + ";Time;Gate", 1000, 0, 10000,
110  192, 0,
111  192);
112  }
113  } else {
114  hOccModAfterInjLER[avxdid] = new TH1F("PXDOccInjLER_" + bufful,
115  "PXDOccModInjLER " + buff + "/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
116  hOccModAfterInjHER[avxdid] = new TH1F("PXDOccInjHER_" + bufful,
117  "PXDOccModInjHER " + buff + "/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
118  if (m_createMaxHist) {
119  hMaxOccModAfterInjLER[avxdid] = new TH1F("PXDMaxOccInjLER_" + bufful,
120  "PXDMaxOccModInjLER " + buff + "/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
121  hMaxOccModAfterInjHER[avxdid] = new TH1F("PXDMaxOccInjHER_" + bufful,
122  "PXDMaxOccModInjHER " + buff + "/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
123  }
124  if (m_createGateHist) {
125  hOccModAfterInjLERGate[avxdid] = new TH2F("PXDOccInjLERGate_" + bufful, "PXDOccInjLERGate " + buff + ";Time;Gate", 1000, 0, 10000,
126  192, 0,
127  192);
128  hOccModAfterInjHERGate[avxdid] = new TH2F("PXDOccInjHERGate_" + bufful, "PXDOccInjHERGate " + buff + ";Time;Gate", 1000, 0, 10000,
129  192, 0,
130  192);
131  }
132  }
133  }
134  }
135 
136 // hTrigAfterInjLER = new TH2F("TrigAfterInjLER",
137 // "Triggers for LER veto tuning;Time since last injection in #mus;Time within beam cycle in #mus", 500, 0, 30000, 100, 0,
138 // 5120 / 508.);
139 // hTrigAfterInjHER = new TH2F("TrigAfterInjHER",
140 // "Triggers for HER veto tuning;Time since last injection in #mus;Time within beam cycle in #mus", 500, 0, 30000, 100, 0,
141 // 5120 / 508.);
142 
143  hTriggersAfterTrigger = new TH1I("PXDTriggersAfterLast",
144  "PXD Trigger after Last Trigger;Time diff in #mus;Count/Time (0.5 #mus bins)", 100000, 0, 50000);
145  hTriggersPerBunch = new TH1I("PXDTriggerBunch", "PXD Trigger per Bunch;Bunch/4;Triggers", 1280, 0, 1280);
146 
147 
148  // cd back to root directory
149  oldDir->cd();
150 }
151 
153 {
154  REG_HISTOGRAM
155  m_EventLevelTriggerTimeInfo.isRequired();
156 
157  if (m_useClusters) {
159  } else {
160  m_storeRawHits.isRequired(m_PXDRawHitsName);
161  }
162 }
163 
165 {
166  // Do not assume that everthing is non-zero, e.g. Max might be nullptr
167  if (hOccAfterInjLER) hOccAfterInjLER->Reset();
168  if (hOccAfterInjHER) hOccAfterInjHER->Reset();
169  if (hEOccAfterInjLER) hEOccAfterInjLER->Reset();
170  if (hEOccAfterInjHER) hEOccAfterInjHER->Reset();
175  for (auto& a : hOccModAfterInjLER) if (a.second) a.second->Reset();
176  for (auto& a : hOccModAfterInjHER) if (a.second) a.second->Reset();
177  for (auto& a : hMaxOccModAfterInjLER) if (a.second) a.second->Reset();
178  for (auto& a : hMaxOccModAfterInjHER) if (a.second) a.second->Reset();
179  for (auto& a : hOccModAfterInjLERGate) if (a.second) a.second->Reset();
180  for (auto& a : hOccModAfterInjHERGate) if (a.second) a.second->Reset();
181 // hTrigAfterInjLER->Reset();
182 // hTrigAfterInjHER->Reset();
183  hTriggersAfterTrigger->Reset();
184  hTriggersPerBunch->Reset();
185 }
186 
188 {
189  // And check if the stored data is valid
190  if (m_EventLevelTriggerTimeInfo.isValid() and m_EventLevelTriggerTimeInfo->isValid()) {
191  // get last injection time
192  hTriggersAfterTrigger->Fill(m_EventLevelTriggerTimeInfo->getTimeSincePrevTrigger() / 127.);
193  // hTriggersAfterTrigger->Fill(m_EventLevelTriggerTimeInfo->getTimeSincePrevTrigger() / 64.);
194  hTriggersPerBunch->Fill(m_EventLevelTriggerTimeInfo->getBunchNumber());
195 
196  // check time overflow, too long ago
197  if (m_EventLevelTriggerTimeInfo->hasInjection()) {
198  // count raw pixel hits or clusters per module, only if necessary
199  unsigned int all = 0;
200  std::map <VxdID, int> freq;// count the number of RawHits per sensor
201  if (m_useClusters) {
202  for (auto& p : m_storeClusters) {
203  freq[p.getSensorID()]++;
204  all++;
205  }
206  } else {
207  for (auto& p : m_storeRawHits) {
208  freq[p.getSensorID()]++;
209  all++;
210  }
211  }
212  float difference = m_EventLevelTriggerTimeInfo->getTimeSinceLastInjection() / 127.; // 127MHz clock ticks to us, inexact rounding
213  if (m_EventLevelTriggerTimeInfo->isHER()) {
214  hOccAfterInjHER->Fill(difference, all);
215  hEOccAfterInjHER->Fill(difference);
216  // hTrigAfterInjHER->Fill(difference, difference - int(difference / (5120 / 508.)) * (5120 / 508.));
217  if (m_createMaxHist) {
218  auto bin = hMaxOccAfterInjHER->FindBin(difference);
219  auto value = hMaxOccAfterInjHER->GetBinContent(bin);
220  if (all > value) hMaxOccAfterInjHER->SetBinContent(bin, all);
221  }
222  for (auto& a : hOccModAfterInjHER) {
223  if (a.second) a.second->Fill(difference, freq[a.first]);
224  }
225  if (m_createMaxHist) {
226  for (auto& a : hMaxOccModAfterInjHER) {
227  if (a.second) {
228  auto bin = a.second->FindBin(difference);
229  auto value = a.second->GetBinContent(bin);
230  if (freq[a.first] > value) a.second->SetBinContent(bin, freq[a.first]);
231  }
232  }
233  }
234  if (hOccAfterInjHERGate) {
235  if (m_useClusters) {
236  // Cluster does not contain VCellID, need to change histogramm completely
237  // -> doesnt work with clusters!
238  // for (auto& p : m_storeClusters) {
239  // hOccAfterInjHERGate->Fill(difference, p.getVCellID() / 4);
240  // }
241  } else {
242  for (auto& p : m_storeRawHits) {
243  hOccAfterInjHERGate->Fill(difference, p.getRow() / 4);
244  hOccModAfterInjHERGate[p.getSensorID()]->Fill(difference, p.getRow() / 4);
245  }
246  }
247  }
248  } else {
249  hOccAfterInjLER->Fill(difference, all);
250  hEOccAfterInjLER->Fill(difference);
251  // hTrigAfterInjLER->Fill(difference, difference - int(difference / (5120 / 508.)) * (5120 / 508.));
252  if (m_createMaxHist) {
253  auto bin = hMaxOccAfterInjLER->FindBin(difference);
254  auto value = hMaxOccAfterInjLER->GetBinContent(bin);
255  if (all > value) hMaxOccAfterInjLER->SetBinContent(bin, all);
256  }
257  for (auto& a : hOccModAfterInjLER) {
258  if (a.second) a.second->Fill(difference, freq[a.first]);
259  }
260  if (m_createMaxHist) {
261  for (auto& a : hMaxOccModAfterInjLER) {
262  if (a.second) {
263  auto bin = a.second->FindBin(difference);
264  auto value = a.second->GetBinContent(bin);
265  if (freq[a.first] > value) a.second->SetBinContent(bin, freq[a.first]);
266  }
267  }
268 
269  }
270  if (hOccAfterInjLERGate) {
271  if (m_useClusters) {
272  // Cluster does not contain VCellID, need to change histogramm completely
273  // -> doesnt work with clusters!
274  // for (auto& p : m_storeClusters) {
275  // hOccAfterInjLERGate->Fill(difference, p.getVCellID() / 4);
276  // }
277  } else {
278  for (auto& p : m_storeRawHits) {
279  hOccAfterInjLERGate->Fill(difference, p.getRow() / 4);
280  hOccModAfterInjLERGate[p.getSensorID()]->Fill(difference, p.getRow() / 4);
281  }
282  }
283  }
284  }
285  }
286  }
287 }
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
std::map< VxdID, TH1F * > hOccModAfterInjLER
Histogram Occupancy after LER injection.
TH1F * hMaxOccAfterInjHER
Histogram Max Occupancy after HER injection.
TH2F * hOccAfterInjLERGate
Occupancy after LER injection per Gate.
std::string m_PXDRawHitsName
The name of the StoreArray of PXDRawHits.
TH1F * hMaxOccAfterInjLER
Histogram Max Occupancy after LER injection.
void initialize() override final
initialize function
TH1I * hTriggersAfterTrigger
Histogram for Nr Entries (=Triggers after Last Trigger.
StoreArray< PXDCluster > m_storeClusters
Input array for PXD Clusters.
std::map< VxdID, TH1F * > hMaxOccModAfterInjHER
Histogram Max Occupancy after HER injection.
std::map< VxdID, TH1F * > hMaxOccModAfterInjLER
Histogram Max Occupancy after LER injection.
void defineHisto() override final
defineHisto function
std::map< VxdID, TH2F * > hOccModAfterInjLERGate
Occupancy after LER injection per Gate per Module.
bool m_eachModule
create a histo per module
bool m_useClusters
use PXDClusters instead of Raw Hits
void event() override final
event function
bool m_createMaxHist
create max hits histogram, not multi processing save!!
std::string m_PXDClustersName
The name of the StoreArray of PXDClusters.
std::map< VxdID, TH2F * > hOccModAfterInjHERGate
Occupancy after HER injection per Gate per Module.
std::string m_histogramDirectoryName
Name of the histogram directory in ROOT file.
TH1F * hOccAfterInjLER
Histogram Occupancy after LER injection.
void beginRun() override final
beginRun function
VXD::GeoCache & m_vxdGeometry
the VXD geometry
TH1F * hOccAfterInjHER
Histogram Occupancy after HER injection.
TH2F * hOccAfterInjHERGate
Occupancy after HER injection per Gate.
TH1I * hEOccAfterInjLER
Histogram for Nr Entries (=Triggrs) for normalization after LER injection.
StoreObjPtr< EventLevelTriggerTimeInfo > m_EventLevelTriggerTimeInfo
Object for TTD mdst object.
StoreArray< PXDRawHit > m_storeRawHits
Input array for PXD Raw Hits.
TH1I * hTriggersPerBunch
Histogram forTrigger per Bunch
std::map< VxdID, TH1F * > hOccModAfterInjHER
Histogram Occupancy after HER injection.
bool m_offlineStudy
create histos with much finer binning and larger range
bool m_createGateHist
create per gate hits 2d histogram
TH1I * hEOccAfterInjHER
Histogram for Nr Entries (=Triggrs) for normalization after HER injection.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
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
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
Namespace to encapsulate code needed for simulation and reconstrucion of the PXD.
Namespace to provide code needed by both Vertex Detectors, PXD and SVD, and also testbeam telescopes.
Definition: GeoCache.h:34
Abstract base class for different kinds of events.