Belle II Software development
PXDDAQDQMModule.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/PXDDAQDQMModule.h>
10
11#include <TDirectory.h>
12#include <TAxis.h>
13#include <boost/format.hpp>
14#include <string>
15
16using namespace std;
17using namespace Belle2;
18using namespace Belle2::PXD;
19using namespace Belle2::PXD::PXDError;
20using namespace Belle2::VXD;
21using boost::format;
22
23//-----------------------------------------------------------------
24// Register the Module
25//-----------------------------------------------------------------
26REG_MODULE(PXDDAQDQM);
27
28//-----------------------------------------------------------------
29// Implementation
30//-----------------------------------------------------------------
31
32PXDDAQDQMModule::PXDDAQDQMModule() : HistoModule(), m_vxdGeometry(VXD::GeoCache::getInstance())
33{
34 //Set module properties
35 setDescription("Monitor DAQ errors");
37 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
38 std::string("pxdDAQ"));
39}
40
42{
43 TDirectory* oldDir = gDirectory;
44 if (m_histogramDirectoryName != "") {
45 oldDir->mkdir(m_histogramDirectoryName.c_str());// do not rely on return value, might be ZERO
46 oldDir->cd(m_histogramDirectoryName.c_str());
47 }
48
49 hDAQErrorEvent = new TH1D("PXDDAQError", "PXDDAQError/Event;;Count", ONSEN_USED_TYPE_ERR, 0, ONSEN_USED_TYPE_ERR);
50 hDAQErrorDHC = new TH2D("PXDDAQDHCError", "PXDDAQError/DHC;DHC ID;", 16, 0, 16, ONSEN_USED_TYPE_ERR, 0, ONSEN_USED_TYPE_ERR);
51 hDAQErrorDHE = new TH2D("PXDDAQDHEError", "PXDDAQError/DHE;DHE ID;", 64, 0, 64, ONSEN_USED_TYPE_ERR, 0, ONSEN_USED_TYPE_ERR);
52 hDAQUseableModule = new TH1D("PXDDAQUseableModule", "PXDDAQUseableModule/DHE;DHE ID;", 64, 0, 64);
53 hDAQNotUseableModule = new TH1D("PXDDAQNotUseableModule", "PXDDAQNotUseableModule/DHE;DHE ID;", 64, 0, 64);
54 hDAQDHPDataMissing = new TH1D("PXDDAQDHPDataMissing", "PXDDAQDHPDataMissing/DHE*DHP;DHE ID;", 64 * 4, 0, 64);
55 hDAQEndErrorDHC = new TH2D("PXDDAQDHCEndError", "PXDDAQEndError/DHC;DHC ID;", 16, 0, 16, 32, 0, 32);
56 hDAQEndErrorDHE = new TH2D("PXDDAQDHEEndError", "PXDDAQEndError/DHE;DHE ID;", 64, 0, 64, 4 * 2 * 8, 0, 4 * 2 * 8);
57
58 // histograms might get unreadable, but, if necessary, you can zoom in anyways.
59 // we could use full alphanumeric histograms, but then, the labels would change (in the worst case) depending on observed errors
60 // and ... the histogram would contain NO labels if there is NO error ... confusing.
61 // ... an we would have to use alphanumeric X axis (DHE ID, DHC ID), too)
62 for (int i = 0; i < ONSEN_USED_TYPE_ERR; i++) {
63 const char* label = getPXDBitErrorName(i).c_str();
64 hDAQErrorEvent->GetXaxis()->SetBinLabel(i + 1, label);
65 hDAQErrorDHE->GetYaxis()->SetBinLabel(i + 1, label);
66 hDAQErrorDHC->GetYaxis()->SetBinLabel(i + 1, label);
67 }
68
69 hDAQErrorEvent->LabelsOption("v"); // rotate the labels.
70
71 std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
72 for (VxdID& avxdid : sensors) {
74 if (info.getType() != VXD::SensorInfoBase::PXD) continue;
75 //Only interested in PXD sensors
76
77 TString buff = (std::string)avxdid;
78 TString bufful = buff;
79 buff.ReplaceAll(".", "_");
80
81// string s = str(format("DHE %d:%d:%d (DHH ID %02Xh)") % num1 % num2 % num3 % i);
82// string s2 = str(format("_%d.%d.%d") % num1 % num2 % num3);
83
84 hDAQDHETriggerGate[avxdid] = new TH1D("PXDDAQDHETriggerGate_" + bufful,
85 "TriggerGate DHE " + buff + "; Trigger Gate; Counts", 192, 0, 192);
86 hDAQDHEReduction[avxdid] = new TH1D("PXDDAQDHEDataReduction_" + bufful, "Data Reduction DHE " + buff + "; Raw/Red; Counts", 200, 0,
87 40);// If max changed, check overflow copy below
88 hDAQCM[avxdid] = new TH2D("PXDDAQCM_" + bufful, "Common Mode on DHE " + buff + "; Gate+Chip*192; Common Mode", 192 * 4, 0, 192 * 4,
89 64, 0, 64);
90 hDAQCM2[avxdid] = new TH1D("PXDDAQCM2_" + bufful, "Common Mode on DHE " + buff + "; Common Mode", 64, 0, 64);
91 }
92 for (int i = 0; i < 16; i++) {
93 hDAQDHCReduction[i] = new TH1D(("PXDDAQDHCDataReduction_" + str(format("%d") % i)).c_str(),
94 ("Data Reduction DHC " + str(format(" %d") % i) + "; Raw/Red; Counts").c_str(), 200, 0,
95 40); // If max changed, check overflow copy below
96 }
97// hDAQErrorEvent->LabelsDeflate("X");
98// hDAQErrorEvent->LabelsOption("v");
99// hDAQErrorEvent->SetStats(0);
100 hEODBAfterInjLER = new TH1I("PXDEODBInjLER", "PXDEODBInjLER/Time;Time in #mus;Events/Time (5 #mus bins)", 4000, 0, 20000);
101 hEODBAfterInjHER = new TH1I("PXDEODBInjHER", "PXDEODBInjHER/Time;Time in #mus;Events/Time (5 #mus bins)", 4000, 0, 20000);
102 hCM63AfterInjLER = new TH1I("PXDCM63InjLER", "PXDCM63InjLER/Time;Time in #mus;Events/Time (5 #mus bins)", 4000, 0, 20000);
103 hCM63AfterInjHER = new TH1I("PXDCM63InjHER", "PXDCM63InjHER/Time;Time in #mus;Events/Time (5 #mus bins)", 4000, 0, 20000);
104 hTruncAfterInjLER = new TH1I("PXDTruncInjLER", "PXDTruncInjLER/Time;Time in #mus;Events/Time (5 #mus bins)", 4000, 0, 20000);
105 hTruncAfterInjHER = new TH1I("PXDTruncInjHER", "PXDTruncInjHER/Time;Time in #mus;Events/Time (5 #mus bins)", 4000, 0, 20000);
106 hMissAfterInjLER = new TH1I("PXDMissInjLER", "PXDMissInjLER/Time;Time in #mus;Events/Time (5 #mus bins)", 4000, 0, 20000);
107 hMissAfterInjHER = new TH1I("PXDMissInjHER", "PXDMissInjHER/Time;Time in #mus;Events/Time (5 #mus bins)", 4000, 0, 20000);
108 hEODBTrgDiff = new TH1I("PXDEODBTrgDiff", "PXDEODBTrgDiff/DiffTime;DiffTime in #mus;Events/Time (1 #mus bins)", 2000, 0, 2000);
109 hCM63TrgDiff = new TH1I("PXDCM63TrgDiff", "PXDCM63TrgDiff/DiffTime;DiffTime in #mus;Events/Time (1 #mus bins)", 2000, 0, 2000);
110 hTruncTrgDiff = new TH1I("PXDTruncTrgDiff", "PXDTruncTrgDiff/DiffTime;DiffTime in #mus;Events/Time (1 #mus bins)", 2000, 0, 2000);
111 hMissTrgDiff = new TH1I("PXDMissTrgDiff", "PXDMissTrgDiff/DiffTime;DiffTime in #mus;Events/Time (1 #mus bins)", 2000, 0, 2000);
112
113 hDAQStat = new TH1D("PXDDAQStat", "PXDDAQStat", 20, 0, 20);
114 auto xa = hDAQStat->GetXaxis();
115 if (xa) {
116 // underflow: number of events -> for normalize
117 xa->SetBinLabel(0 + 1, "EODB/HLT rej"); // event of doom or rejected
118 xa->SetBinLabel(1 + 1, "Trunc 8%");
119 xa->SetBinLabel(2 + 1, "HER Trunc");
120 xa->SetBinLabel(3 + 1, "LER Trunc");
121 xa->SetBinLabel(4 + 1, "CM63");
122 xa->SetBinLabel(5 + 1, "HER CM63");
123 xa->SetBinLabel(6 + 1, "LER CM63");
124 xa->SetBinLabel(7 + 1, "HER CM63>1ms");
125 xa->SetBinLabel(8 + 1, "LER CM63>1ms");
126 xa->SetBinLabel(9 + 1, "HER Trunc>1ms");
127 xa->SetBinLabel(10 + 1, "LER Trunc>1ms");
128 xa->SetBinLabel(11 + 1, "MissFrame");
129 xa->SetBinLabel(12 + 1, "Timeout");
130 xa->SetBinLabel(13 + 1, "Link Down");
131 xa->SetBinLabel(14 + 1, "Mismatch");
132 xa->SetBinLabel(15 + 1, "HER MissFrame");
133 xa->SetBinLabel(16 + 1, "LER MissFrame");
134 xa->SetBinLabel(17 + 1, "HER MissFrame>1ms");
135 xa->SetBinLabel(18 + 1, "LER MissFrame>1ms");
136 xa->SetBinLabel(19 + 1, "unused");
137 }
138 // cd back to root directory
139 oldDir->cd();
140}
141
143{
144 REG_HISTOGRAM
145 m_storeDAQEvtStats.isRequired();
146 m_rawSVD.isOptional();
147 m_EventLevelTriggerTimeInfo.isRequired();
148
149}
150
152{
153 hDAQErrorEvent->Reset();
154 hDAQErrorDHC->Reset();
155 hDAQErrorDHE->Reset();
156 hDAQUseableModule->Reset();
157 hDAQNotUseableModule->Reset();
158 hDAQDHPDataMissing->Reset();
159 hDAQEndErrorDHC->Reset();
160 hDAQEndErrorDHE->Reset();
161 for (auto& it : hDAQDHETriggerGate) if (it.second) it.second->Reset();
162 for (auto& it : hDAQDHCReduction) if (it.second) it.second->Reset();
163 for (auto& it : hDAQDHEReduction) if (it.second) it.second->Reset();
164 for (auto& it : hDAQCM) if (it.second) it.second->Reset();
165 for (auto& it : hDAQCM2) if (it.second) it.second->Reset();
172 hDAQStat->Reset();
173}
174
176{
177 hDAQErrorEvent->Fill(-1);// Event counter
178 hDAQStat->Fill(-1); // to normalize to the number of events
179 hDAQDHPDataMissing->Fill(-1); // to normalize to the number of events
180 hDAQErrorDHC->Fill(-1, -1); // to normalize to the number of events
181 hDAQErrorDHE->Fill(-1, -1); // to normalize to the number of events
182 for (auto& it : hDAQCM2) if (it.second) it.second->Fill(-1); // to normalize to the number of events
183 for (auto& it : hDAQCM) if (it.second) it.second->Fill(-1, -1); // to normalize to the number of events
188
189 bool eodbFlag = m_rawSVD.getEntries() == 0;
190
191 bool truncFlag = false; // flag events which are DHE truncated
192 bool nolinkFlag = false; // flag events which are DHE truncated
193 bool missingFlag = false; // flag events where frame is missing
194 bool timeoutFlag = false; // flag events where frame timeout
195 bool mismatchFlag = false; // flag events where trig mismatched
196 bool cm63Flag = false; // flag event which are CM63 truncated
197
198 B2DEBUG(20, "Iterate PXD DAQ Status");
199 auto evt = *m_storeDAQEvtStats;
200 auto evt_emask = evt.getErrorMask();
201 for (int i = 0; i < ONSEN_MAX_TYPE_ERR; i++) {
202 if (evt_emask[i]) hDAQErrorEvent->Fill(getPXDBitErrorName(i).c_str(), 1);
203 }
204 B2DEBUG(20, "Iterate PXD Packets, Err " << evt_emask);
205 for (auto& pkt : evt) {
206 B2DEBUG(20, "Iterate PXD DHC in Pkt " << pkt.getPktIndex());
207 for (auto& dhc : pkt) {
208 hDAQErrorDHC->Fill(dhc.getDHCID(), -1);// normalize
209 auto dhc_emask = dhc.getErrorMask();
210 for (int i = 0; i < ONSEN_MAX_TYPE_ERR; i++) {
211 if (dhc_emask[i]) hDAQErrorDHC->Fill(dhc.getDHCID(), i);
212 }
213 unsigned int cmask = dhc.getEndErrorInfo();
214 for (int i = 0; i < 32; i++) {
215 unsigned int mask = (1u << i);
216 if ((cmask & mask) == mask) hDAQEndErrorDHC->Fill(dhc.getDHCID(), i);
217 }
218 if (hDAQDHCReduction[dhc.getDHCID()]) {
219 float red = dhc.getRedCnt() ? float(dhc.getRawCnt()) / dhc.getRedCnt() : 0.;
220 B2DEBUG(98, "==DHC " << dhc.getDHCID() << "(Raw)" << dhc.getRawCnt() << " / (Red)" << dhc.getRedCnt() << " = " << red);
221 if (red >= 40.) red = 39.999999999; // Bad, bad workaround. but we want to see the overflows
222 hDAQDHCReduction[dhc.getDHCID()]->Fill(red);
223 }
224 B2DEBUG(20, "Iterate PXD DHE in DHC " << dhc.getDHCID() << " , Err " << dhc_emask);
225 for (auto& dhe : dhc) {
226 hDAQErrorDHE->Fill(dhe.getDHEID(), -1);// normalize
227 auto dhe_emask = dhe.getErrorMask();
228 B2DEBUG(20, "DHE " << dhe.getDHEID() << " , Err " << dhe_emask);
229 for (int i = 0; i < ONSEN_MAX_TYPE_ERR; i++) {
230 if (dhe_emask[i]) hDAQErrorDHE->Fill(dhe.getDHEID(), i);
231 }
232 if (dhe.isUsable()) {
233 hDAQUseableModule->Fill(dhe.getDHEID());
234 } else {
235 hDAQNotUseableModule->Fill(dhe.getDHEID());
236 }
237 for (int i = 0; i < 4; i++) {
238 if ((dhe.getDHPFoundMask() & (1 << i)) == 0) hDAQDHPDataMissing->Fill(dhe.getDHEID() + i * 0.25);
239 }
240 for (auto& dhp : dhe) {
241 truncFlag |= dhp.getTruncated(); // new firmware workaround flag
242 }
243 unsigned int emask = dhe.getEndErrorInfo();
244 // TODO differentiate between link-lost and truncation
245 for (int i = 0; i < 4 * 2; i++) {
246 auto sm = (emask >> i * 4) & 0xF;
247 if (sm >= 8) sm = 7; // clip unknow to 7, as value >6 undefined for now
248 if (sm > 0) hDAQEndErrorDHE->Fill(dhe.getDHEID(), i * 8 + sm); // we dont want to fill noerror=0
249 missingFlag |= sm == 0x1; // missing
250 timeoutFlag |= sm == 0x2; // timeout
251 nolinkFlag |= sm == 0x3; // link down
252 // 4 is DHP masked
253 mismatchFlag |= sm == 0x5; // start/end mismatch
254 truncFlag |= sm == 0x6; // trunc because of size
255 }
256
257 if (hDAQDHETriggerGate[dhe.getSensorID()]) hDAQDHETriggerGate[dhe.getSensorID()]->Fill(dhe.getTriggerGate());
258 if (hDAQDHEReduction[dhe.getSensorID()]) {
259 float red = dhe.getRedCnt() ? float(dhe.getRawCnt()) / dhe.getRedCnt() : 0.;
260 B2DEBUG(98, "==DHE " << dhe.getSensorID() << "(Raw)" << dhe.getRawCnt() << " / (Red)" << dhe.getRedCnt() << " = " << red);
261 if (red >= 40.) red = 39.999999999; // Bad, bad workaround. but we want to see the overflows
262 hDAQDHEReduction[dhe.getSensorID()]->Fill(red);
263 }
264 for (auto cm = dhe.cm_begin(); cm < dhe.cm_end(); ++cm) {
265 // uint8_t, uint16_t, uint8_t ; tuple of Chip ID (2 bit), Row (10 bit), Common Mode (6 bit)
266 if (hDAQCM[dhe.getSensorID()]) hDAQCM[dhe.getSensorID()]->Fill(std::get<0>(*cm) * 192 + std::get<1>(*cm) / 4, std::get<2>(*cm));
267 if (hDAQCM2[dhe.getSensorID()]) hDAQCM2[dhe.getSensorID()]->Fill(std::get<2>(*cm));
268 //if (hDAQCM2[dhe.getSensorID()] && std::get<1>(*cm) < 768 - 2) {
269 // Deactivated again as the clean-up effect was insignificant for beam data
270 // ignore the always-on rows (second to last: 766 and last 767)
271 // row in CM array is already in readout-direction and not V-ID
272 // hDAQCM2[dhe.getSensorID()]->Fill(std::get<2>(*cm));
273 //}
274 cm63Flag |= 63 == std::get<2>(*cm);
275 }
276 }
277 }
278 }
279
280
281
282 // And check if the stored data is valid
283 if (m_EventLevelTriggerTimeInfo.isValid() and m_EventLevelTriggerTimeInfo->isValid()) {
284
285 double lasttrig = m_EventLevelTriggerTimeInfo->getTimeSincePrevTrigger() / 127.; // 127MHz clock ticks to us, inexact rounding
286 if (eodbFlag && hEODBTrgDiff) hEODBTrgDiff->Fill(lasttrig);
287 if (cm63Flag && hCM63TrgDiff) hCM63TrgDiff->Fill(lasttrig);
288 if (truncFlag && hTruncTrgDiff) hTruncTrgDiff->Fill(lasttrig);
289 if (missingFlag && hMissTrgDiff) hMissTrgDiff->Fill(lasttrig);
290
291 // check time overflow, too long ago
292 if (m_EventLevelTriggerTimeInfo->hasInjection()) {
293 // get last injection time
294 double difference = m_EventLevelTriggerTimeInfo->getTimeSinceLastInjection() / 127.; // 127MHz clock ticks to us, inexact rounding
295 if (m_EventLevelTriggerTimeInfo->isHER()) {
296 if (eodbFlag) {
297 if (hEODBAfterInjHER) hEODBAfterInjHER->Fill(difference);
298 }
299 if (cm63Flag) {
300 hDAQStat->Fill(5); // sum CM63 after HER
301 if (difference > 1000) hDAQStat->Fill(7); // sum CM63 after HER, but outside injections, 1ms
302 if (hCM63AfterInjHER) hCM63AfterInjHER->Fill(difference);
303 }
304 if (truncFlag) {
305 hDAQStat->Fill(2); // sum truncs after HER
306 if (difference > 1000) hDAQStat->Fill(9); // sum truncs after HER, but outside injections, 1ms
307 if (hTruncAfterInjHER) hTruncAfterInjHER->Fill(difference);
308 }
309 if (missingFlag) {
310 hDAQStat->Fill(15); // sum missframe after HER
311 if (difference > 1000) hDAQStat->Fill(17); // sum missframe after HER, but outside injections, 1ms
312 if (hMissAfterInjHER) hMissAfterInjHER->Fill(difference);
313 }
314 } else {
315 if (eodbFlag) {
316 if (hEODBAfterInjLER) hEODBAfterInjLER->Fill(difference);
317 }
318 if (cm63Flag) {
319 hDAQStat->Fill(6); // sum CM63 after LER
320 if (difference > 1000) hDAQStat->Fill(8); // sum CM63 after LER, but outside injections, 1ms
321 if (hCM63AfterInjLER) hCM63AfterInjLER->Fill(difference);
322 }
323 if (truncFlag) {
324 hDAQStat->Fill(3); // sum truncs after LER
325 if (difference > 1000) hDAQStat->Fill(10); // sum truncs after LER, but outside injections, 1ms
326 if (hTruncAfterInjLER) hTruncAfterInjLER->Fill(difference);
327 }
328 if (missingFlag) {
329 hDAQStat->Fill(16); // sum missframe after LER
330 if (difference > 1000) hDAQStat->Fill(18); // sum missframe after LER, but outside injections, 1ms
331 if (hMissAfterInjLER) hMissAfterInjLER->Fill(difference);
332 }
333 }
334 }
335 }
336
337// make some nice statistics
338 if (truncFlag) hDAQStat->Fill(1);
339 if (cm63Flag) hDAQStat->Fill(4);
340 if (missingFlag) hDAQStat->Fill(11);
341 if (timeoutFlag) hDAQStat->Fill(12);
342 if (nolinkFlag) hDAQStat->Fill(13);
343 if (mismatchFlag) hDAQStat->Fill(14);
344
345// Check Event-of-doom-busted or otherwise HLT rejected events
346 if (eodbFlag) hDAQStat->Fill(0);
347}
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
TH1I * hCM63AfterInjLER
Histogram of CM63 after LER injection.
void initialize() override final
Initialize.
PXDDAQDQMModule()
Constructor defining the parameters.
std::map< VxdID, TH1D * > hDAQDHETriggerGate
DHE Trigger Gate ("start Row")
TH1D * hDAQNotUseableModule
Count Useable/unuseable decision.
StoreObjPtr< PXDDAQStatus > m_storeDAQEvtStats
Input array for DAQ Status.
TH1I * hTruncAfterInjHER
Histogram Truncation after HER injection.
StoreArray< RawSVD > m_rawSVD
Input array for SVD/x-check HLT EODB .
TH1D * hDAQErrorEvent
Remark: Because of DHH load balancing and sub event building, the very same DHE and DHC can show up i...
TH1I * hEODBTrgDiff
Histogram of EODB after last trigger.
std::map< VxdID, TH2D * > hDAQCM
Common Mode per DHE to gate and DHP level.
TH1D * hDAQUseableModule
Count Useable/unuseable decision.
void defineHisto() override final
Define histograms.
TH1I * hMissTrgDiff
Histogram MissFrame after last trigger.
TH1I * hMissAfterInjHER
Histogram MissFrame after HER injection.
void event() override final
Event.
TH2D * hDAQErrorDHC
individual DHC errors
TH2D * hDAQEndErrorDHE
individual DHE END errors
std::string m_histogramDirectoryName
Name of the histogram directory in ROOT file.
TH1I * hCM63TrgDiff
Histogram of CM63 after last trigger.
TH1I * hCM63AfterInjHER
Histogram of CM63 after HER injection.
TH1I * hTruncTrgDiff
Histogram Truncation after last trigger.
std::map< VxdID, TH1D * > hDAQDHEReduction
DHE data reduction
TH1I * hEODBAfterInjHER
Histogram of EODB after HER injection.
TH2D * hDAQEndErrorDHC
individual DHC END errors
TH2D * hDAQErrorDHE
individual DHE errors
void beginRun() override final
Begin run.
std::map< VxdID, TH1D * > hDAQCM2
Common Mode per DHE to gate and DHP level.
VXD::GeoCache & m_vxdGeometry
the geometry
TH1I * hMissAfterInjLER
Histogram MissFrame after LER injection.
TH1D * hDAQDHPDataMissing
Count Missing DHP data.
std::map< int, TH1D * > hDAQDHCReduction
DHC data reduction
TH1D * hDAQStat
Histogram for Truncation etc Stats.
TH1I * hTruncAfterInjLER
Histogram Truncation after LER injection.
StoreObjPtr< EventLevelTriggerTimeInfo > m_EventLevelTriggerTimeInfo
Object for TTD mdst object.
TH1I * hEODBAfterInjLER
Histogram of EODB after LER 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: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.
STL namespace.