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