Belle II Software  release-06-02-00
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");
30  setPropertyFlags(c_ParallelProcessingCertified);
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 
43 void PXDInjectionDQMModule::defineHisto()
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) {
84  VXD::SensorInfoBase info = m_vxdGeometry.getSensorInfo(avxdid);
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 
152 void PXDInjectionDQMModule::initialize()
153 {
154  REG_HISTOGRAM
155  m_rawTTD.isOptional();
156  if (m_useClusters) {
157  m_storeClusters.isRequired(m_PXDClustersName);
158  } else {
159  m_storeRawHits.isRequired(m_PXDRawHitsName);
160  }
161 }
162 
163 void PXDInjectionDQMModule::beginRun()
164 {
165  // Do not assume that everthing is non-zero, e.g. Max might be nullptr
166  if (hOccAfterInjLER) hOccAfterInjLER->Reset();
167  if (hOccAfterInjHER) hOccAfterInjHER->Reset();
168  if (hEOccAfterInjLER) hEOccAfterInjLER->Reset();
169  if (hEOccAfterInjHER) hEOccAfterInjHER->Reset();
170  if (hMaxOccAfterInjLER) hMaxOccAfterInjLER->Reset();
171  if (hMaxOccAfterInjHER) hMaxOccAfterInjHER->Reset();
172  if (hOccAfterInjLERGate) hOccAfterInjLERGate->Reset();
173  if (hOccAfterInjHERGate) hOccAfterInjHERGate->Reset();
174  for (auto& a : hOccModAfterInjLER) if (a.second) a.second->Reset();
175  for (auto& a : hOccModAfterInjHER) if (a.second) a.second->Reset();
176  for (auto& a : hMaxOccModAfterInjLER) if (a.second) a.second->Reset();
177  for (auto& a : hMaxOccModAfterInjHER) if (a.second) a.second->Reset();
178  for (auto& a : hOccModAfterInjLERGate) if (a.second) a.second->Reset();
179  for (auto& a : hOccModAfterInjHERGate) if (a.second) a.second->Reset();
180 // hTrigAfterInjLER->Reset();
181 // hTrigAfterInjHER->Reset();
182  hTriggersAfterTrigger->Reset();
183  hTriggersPerBunch->Reset();
184 }
185 
186 void PXDInjectionDQMModule::event()
187 {
188 
189  for (auto& it : m_rawTTD) {
190  B2DEBUG(29, "TTD FTSW : " << hex << it.GetTTUtime(0) << " " << it.GetTTCtime(0) << " EvtNr " << it.GetEveNo(0) << " Type " <<
191  (it.GetTTCtimeTRGType(0) & 0xF) << " TimeSincePrev " << it.GetTimeSincePrevTrigger(0) << " TimeSinceInj " <<
192  it.GetTimeSinceLastInjection(0) << " IsHER " << it.GetIsHER(0) << " Bunch " << it.GetBunchNumber(0));
193 
194  // get last injection time
195  hTriggersAfterTrigger->Fill(it.GetTimeSincePrevTrigger(0) / 127.);
196  hTriggersPerBunch->Fill(it.GetBunchNumber(0));
197 
198  auto difference = it.GetTimeSinceLastInjection(0);
199  // check time overflow, too long ago
200  if (difference != 0x7FFFFFFF) {
201  // count raw pixel hits or clusters per module, only if necessary
202  unsigned int all = 0;
203  std::map <VxdID, int> freq;// count the number of RawHits per sensor
204  if (m_useClusters) {
205  for (auto& p : m_storeClusters) {
206  freq[p.getSensorID()]++;
207  all++;
208  }
209  } else {
210  for (auto& p : m_storeRawHits) {
211  freq[p.getSensorID()]++;
212  all++;
213  }
214  }
215  float diff2 = difference / 127.; // 127MHz clock ticks to us, inexact rounding
216  if (it.GetIsHER(0)) {
217  hOccAfterInjHER->Fill(diff2, all);
218  hEOccAfterInjHER->Fill(diff2);
219 // hTrigAfterInjHER->Fill(diff2, diff2 - int(diff2 / (5120 / 508.)) * (5120 / 508.));
220  if (m_createMaxHist) {
221  auto bin = hMaxOccAfterInjHER->FindBin(diff2);
222  auto value = hMaxOccAfterInjHER->GetBinContent(bin);
223  if (all > value) hMaxOccAfterInjHER->SetBinContent(bin, all);
224  }
225  for (auto& a : hOccModAfterInjHER) {
226  if (a.second) a.second->Fill(diff2, freq[a.first]);
227  }
228  if (m_createMaxHist) {
229  for (auto& a : hMaxOccModAfterInjHER) {
230  if (a.second) {
231  auto bin = a.second->FindBin(diff2);
232  auto value = a.second->GetBinContent(bin);
233  if (freq[a.first] > value) a.second->SetBinContent(bin, freq[a.first]);
234  }
235  }
236  }
237  if (hOccAfterInjHERGate) {
238  if (m_useClusters) {
239  // Cluster does not contain VCellID, need to change histogramm completely
240  // -> doesnt work with clusters!
241 // for (auto& p : m_storeClusters) {
242 // hOccAfterInjHERGate->Fill(diff2, p.getVCellID() / 4);
243 // }
244  } else {
245  for (auto& p : m_storeRawHits) {
246  hOccAfterInjHERGate->Fill(diff2, p.getRow() / 4);
247  hOccModAfterInjHERGate[p.getSensorID()]->Fill(diff2, p.getRow() / 4);
248  }
249  }
250  }
251  } else {
252  hOccAfterInjLER->Fill(diff2, all);
253  hEOccAfterInjLER->Fill(diff2);
254 // hTrigAfterInjLER->Fill(diff2, diff2 - int(diff2 / (5120 / 508.)) * (5120 / 508.));
255  if (m_createMaxHist) {
256  auto bin = hMaxOccAfterInjLER->FindBin(diff2);
257  auto value = hMaxOccAfterInjLER->GetBinContent(bin);
258  if (all > value) hMaxOccAfterInjLER->SetBinContent(bin, all);
259  }
260  for (auto& a : hOccModAfterInjLER) {
261  if (a.second) a.second->Fill(diff2, freq[a.first]);
262  }
263  if (m_createMaxHist) {
264  for (auto& a : hMaxOccModAfterInjLER) {
265  if (a.second) {
266  auto bin = a.second->FindBin(diff2);
267  auto value = a.second->GetBinContent(bin);
268  if (freq[a.first] > value) a.second->SetBinContent(bin, freq[a.first]);
269  }
270  }
271 
272  }
273  if (hOccAfterInjLERGate) {
274  if (m_useClusters) {
275  // Cluster does not contain VCellID, need to change histogramm completely
276  // -> doesnt work with clusters!
277 // for (auto& p : m_storeClusters) {
278 // hOccAfterInjLERGate->Fill(diff2, p.getVCellID() / 4);
279 // }
280  } else {
281  for (auto& p : m_storeRawHits) {
282  hOccAfterInjLERGate->Fill(diff2, p.getRow() / 4);
283  hOccModAfterInjLERGate[p.getSensorID()]->Fill(diff2, p.getRow() / 4);
284  }
285  }
286  }
287  }
288  }
289 
290  break;
291  }
292 }
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
The PXD Occupancy after Injection DQM module.
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
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
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.