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