Belle II Software development
KLMDQMModule.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/* Own header. */
10#include <klm/modules/KLMDQM/KLMDQMModule.h>
11
12/* KLM headers. */
13#include <klm/dataobjects/KLMChannelIndex.h>
14#include <klm/dataobjects/KLMDigitRaw.h>
15
16/* ROOT headers. */
17#include <TDirectory.h>
18
19using namespace Belle2;
20
21REG_MODULE(KLMDQM);
22
25 m_DAQInclusion{nullptr},
26 m_TimeRPC{nullptr},
27 m_TimeScintillatorBKLM{nullptr},
28 m_TimeScintillatorEKLM{nullptr},
29 m_PlaneBKLMPhi{nullptr},
30 m_PlaneBKLMZ{nullptr},
31 m_PlaneEKLM{nullptr},
32 m_MaskedChannelsPerSector{nullptr},
33 m_DigitsKLM{nullptr},
34 m_DigitsRPC{nullptr},
35 m_DigitsScintillatorBKLM{nullptr},
36 m_DigitsScintillatorEKLM{nullptr},
37 m_DigitsMultiStripBKLM{nullptr},
38 m_DigitsMultiStripEKLM{nullptr},
39 m_TriggerBitsBKLM{nullptr},
40 m_TriggerBitsEKLM{nullptr},
41 m_DigitsAfterLERInj{nullptr},
42 m_TriggersLERInj{nullptr},
43 m_DigitsAfterHERInj{nullptr},
44 m_TriggersHERInj{nullptr},
45 m_ChannelArrayIndex{&(KLMChannelArrayIndex::Instance())},
46 m_SectorArrayIndex{&(KLMSectorArrayIndex::Instance())},
47 m_ElementNumbers{&(KLMElementNumbers::Instance())},
48 m_eklmElementNumbers{&(EKLMElementNumbers::Instance())}
49{
50 setDescription("KLM data quality monitor.");
52 addParam("histogramDirectoryName", m_HistogramDirectoryName,
53 "Directory for KLM DQM histograms in ROOT file.",
54 std::string("KLM"));
55 addParam("RPCTimeMin", m_RPCTimeMin,
56 "Min time for RPC time histogram.", double(-1223.5));
57 addParam("RPCTimeMax", m_RPCTimeMax,
58 "Max time for RPC time histogram.", double(-199.5));
59 addParam("BKLMScintTimeMin", m_BKLMScintTimeMin,
60 "Min time for BKLM scintillator time histogram.", double(-5300));
61 addParam("BKLMScintTimeMax", m_BKLMScintTimeMax,
62 "Max time for BKLM scintillator time histogram.", double(-4300));
63 addParam("EKLMScintTimeMin", m_EKLMScintTimeMin,
64 "Min time for EKLM scintillator time histogram.", double(-5300));
65 addParam("EKLMScintTimeMax", m_EKLMScintTimeMax,
66 "Max time for EKLM scintillator time histogram.", double(-4300));
67}
68
70{
72 for (KLMChannelIndex& klmSector : klmIndex) {
73 KLMSectorNumber sector = klmSector.getKLMSectorNumber();
74 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sector);
75 if (m_ChannelHits[sectorIndex] != nullptr)
76 delete[] m_ChannelHits[sectorIndex];
77 }
79 for (KLMChannelIndex& klmSection : klmIndex) {
80 KLMSubdetectorNumber subdetector = klmSection.getSubdetector();
81 if (subdetector == KLMElementNumbers::c_EKLM) {
82 KLMSectionNumber section = klmSection.getSection();
83 if (m_Spatial2DHitsEKLM[section - 1] != nullptr)
84 delete[] m_Spatial2DHitsEKLM[section - 1];
85 }
86 }
87}
88
90{
91 TDirectory* oldDirectory, *newDirectory;
92 oldDirectory = gDirectory;
93 newDirectory = oldDirectory->mkdir(m_HistogramDirectoryName.c_str());
94 newDirectory->cd();
95 /* DAQ inclusion. */
96 m_DAQInclusion = new TH1F("daq_inclusion", "Is KLM included in DAQ?", 2, 0.0, 2.0);
97 m_DAQInclusion->GetXaxis()->SetBinLabel(1, "No");
98 m_DAQInclusion->GetXaxis()->SetBinLabel(2, "Yes");
99 /* Time histograms. */
100 m_TimeRPC = new TH1F("time_rpc", "RPC hit time", 128, m_RPCTimeMin, m_RPCTimeMax);
101 m_TimeRPC->GetXaxis()->SetTitle("Time, ns");
103 new TH1F("time_scintillator_bklm", "Scintillator hit time (BKLM)",
105 m_TimeScintillatorBKLM->GetXaxis()->SetTitle("Time, ns");
107 new TH1F("time_scintillator_eklm", "Scintillator hit time (EKLM)",
109 m_TimeScintillatorEKLM->GetXaxis()->SetTitle("Time, ns");
110 /* Number of hits per plane. */
111 m_PlaneBKLMPhi = new TH1F("plane_bklm_phi",
112 "BKLM plane occupancy (#phi readout)",
113 240, 0.5, 240.5);
114 m_PlaneBKLMPhi->GetXaxis()->SetTitle("Layer number");
115 m_PlaneBKLMZ = new TH1F("plane_bklm_z",
116 "BKLM plane occupancy (z readout)",
117 240, 0.5, 240.5);
118 m_PlaneBKLMZ->GetXaxis()->SetTitle("Layer number");
119 m_PlaneEKLM = new TH1F("plane_eklm", "EKLM plane occupancy (both readouts)", 208, 0.5, 208.5);
120 m_PlaneEKLM->GetXaxis()->SetTitle("Plane number");
121 /* Number of hits per channel. */
122 int nChannelHistograms =
127 KLMChannelNumber* firstChannelNumbers =
128 new KLMChannelNumber[nChannelHistograms + 1];
129 int i = 0;
131 for (KLMChannelIndex& klmSector : klmIndex) {
132 KLMChannelIndex klmChannel(klmSector);
134 KLMChannelNumber channel = klmChannel.getKLMChannelNumber();
135 firstChannelNumbers[i] = m_ChannelArrayIndex->getIndex(channel);
136 if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM) {
138 klmChannel.getSection(), klmChannel.getSector(), 8, 0, 1);
139 firstChannelNumbers[i + 1] = m_ChannelArrayIndex->getIndex(channel);
140 i += 2;
141 } else {
142 int layerIncrease = (klmSector.getSection() == 1) ? 4 : 5;
144 klmChannel.getSection(), klmChannel.getSector(),
145 1 + layerIncrease, 1, 1);
146 firstChannelNumbers[i + 1] = m_ChannelArrayIndex->getIndex(channel);
148 klmChannel.getSection(), klmChannel.getSector(),
149 1 + layerIncrease * 2, 1, 1);
150 firstChannelNumbers[i + 2] = m_ChannelArrayIndex->getIndex(channel);
151 i += 3;
152 }
153 }
154 firstChannelNumbers[nChannelHistograms] = m_ChannelArrayIndex->getNElements();
155 i = 0;
157 for (KLMChannelIndex& klmSector : klmIndex) {
158 int nHistograms;
159 if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM)
160 nHistograms = m_ChannelHitHistogramsBKLM;
161 else
162 nHistograms = m_ChannelHitHistogramsEKLM;
163 KLMSectorNumber sector = klmSector.getKLMSectorNumber();
164 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sector);
165 m_ChannelHits[sectorIndex] = new TH1F*[nHistograms];
166 for (int j = 0; j < nHistograms; j++) {
167 std::string name =
168 "strip_hits_subdetector_" + std::to_string(klmSector.getSubdetector()) +
169 "_section_" + std::to_string(klmSector.getSection()) +
170 "_sector_" + std::to_string(klmSector.getSector()) +
171 "_" + std::to_string(j);
172 std::string title = "Sector " + std::to_string(klmSector.getSector()) + " -- " +
173 m_ElementNumbers->getSectorDAQName(klmSector.getSubdetector(), klmSector.getSection(), klmSector.getSector());
174 m_ChannelHits[sectorIndex][j] = new TH1F(
175 name.c_str(), title.c_str(),
176 firstChannelNumbers[i + 1] - firstChannelNumbers[i],
177 firstChannelNumbers[i] - 0.5, firstChannelNumbers[i + 1] - 0.5);
178 i++;
179 }
180 }
181 delete[] firstChannelNumbers;
182 /* Masked channels per sector:
183 * it is defined here, but filled by the analysis module. */
185 m_MaskedChannelsPerSector = new TH1F("masked_channels", "Number of masked channels per sector",
186 totalSectors, -0.5, totalSectors - 0.5);
188 for (KLMChannelIndex& klmSector : klmIndex) {
189 std::string label = m_ElementNumbers->getSectorDAQName(klmSector.getSubdetector(), klmSector.getSection(), klmSector.getSector());
190 KLMSectorNumber sector = klmSector.getKLMSectorNumber();
191 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sector);
192 m_MaskedChannelsPerSector->GetXaxis()->SetBinLabel(sectorIndex + 1, label.c_str());
193 }
194 /* Number of digits. */
195 m_DigitsKLM = new TH1F("digits_klm", "Number of KLM digits",
196 250.0, 0.0, 250.0);
197 m_DigitsKLM->GetXaxis()->SetTitle("Number of digits");
198 m_DigitsRPC = new TH1F("digits_rpc", "Number of RPC digits",
199 250.0, 0.0, 250.0);
200 m_DigitsRPC->GetXaxis()->SetTitle("Number of digits");
201 m_DigitsScintillatorBKLM = new TH1F("digits_scintillator_bklm", "Number of scintillator digits (BKLM)",
202 250.0, 0.0, 250.0);
203 m_DigitsScintillatorBKLM->GetXaxis()->SetTitle("Number of digits");
204 m_DigitsScintillatorEKLM = new TH1F("digits_scintillator_eklm", "Number of scintillator digits (EKLM)",
205 250.0, 0.0, 250.0);
206 m_DigitsScintillatorEKLM->GetXaxis()->SetTitle("Number of digits");
207 m_DigitsMultiStripBKLM = new TH1F("digits_multi_bklm", "Number of multi-strip digits (BKLM)",
208 50.0, 0.0, 50.0);
209 m_DigitsMultiStripBKLM->GetXaxis()->SetTitle("Number of multi-strip digits");
210 m_DigitsMultiStripEKLM = new TH1F("digits_multi_eklm", "Number of multi-strip digits (EKLM)",
211 50.0, 0.0, 50.0);
212 m_DigitsMultiStripEKLM->GetXaxis()->SetTitle("Number of multi-strip digits");
213 /* Trigger bits. */
214 m_TriggerBitsBKLM = new TH1F("trigger_bits_bklm", "Trigger bits of multi-strip digits (BKLM)",
215 (double)c_0x1, (double)c_0x8, (double)c_0x1 + 1.0);
216 m_TriggerBitsBKLM->GetXaxis()->SetBinLabel(c_0x8, "0x8");
217 m_TriggerBitsBKLM->GetXaxis()->SetBinLabel(c_0x4, "0x4");
218 m_TriggerBitsBKLM->GetXaxis()->SetBinLabel(c_0x2, "0x2");
219 m_TriggerBitsBKLM->GetXaxis()->SetBinLabel(c_0x1, "0x1");
220 m_TriggerBitsEKLM = new TH1F("trigger_bits_eklm", "Trigger bits of multi-strip digits (EKLM)",
221 (double)c_0x1, (double)c_0x8, (double)c_0x1 + 1.0);
222 m_TriggerBitsEKLM->GetXaxis()->SetBinLabel(c_0x8, "0x8");
223 m_TriggerBitsEKLM->GetXaxis()->SetBinLabel(c_0x4, "0x4");
224 m_TriggerBitsEKLM->GetXaxis()->SetBinLabel(c_0x2, "0x2");
225 m_TriggerBitsEKLM->GetXaxis()->SetBinLabel(c_0x1, "0x1");
226 /* Number of digits after injection */
227 /* For the histograms below, we use the same style as for other subdetectors. */
228 m_DigitsAfterLERInj = new TH1F("KLMOccInjLER", "KLM digits after LER injection / Time;Time [#mus];Number of KLM digits / (5 #mus)",
229 4000, 0, 20000);
230 m_TriggersLERInj = new TH1F("KLMTrigInjLER", "Triggers after KER injection / Time;Time [#mus];Number of triggers / (5 #mus)",
231 4000, 0, 20000);
232 m_DigitsAfterHERInj = new TH1F("KLMOccInjHER", "KLM digits after HER injection / Time;Time [#mus];Number of KLM digits / (5 #mus)",
233 4000, 0, 20000);
234 m_TriggersHERInj = new TH1F("KLMTrigInjHER", "Triggers after HER injection / Time;Time [#mus];Number of triggers / (5 #mus)",
235 4000, 0, 20000);
236 /* Spatial distribution of EKLM 2d hits per layer. */
238 for (KLMChannelIndex& klmSection : klmIndex) {
239 KLMSubdetectorNumber subdetector = klmSection.getSubdetector();
240 if (subdetector == KLMElementNumbers::c_EKLM) {
241 KLMSectionNumber section = klmSection.getSection();
242 int maximalLayerNumber = m_eklmElementNumbers->getMaximalDetectorLayerNumber(section);
243 m_Spatial2DHitsEKLM[section - 1] = new TH2F*[maximalLayerNumber];
244 std::string sectionName = (section == EKLMElementNumbers::c_ForwardSection) ? "Forward" : "Backward";
245 for (int j = 1; j <= maximalLayerNumber; ++j) {
246 std::string name = "spatial_2d_hits_subdetector_" + std::to_string(subdetector) +
247 "_section_" + std::to_string(section) +
248 "_layer_" + std::to_string(j);
249 std::string title = "Endcap " + sectionName + " , Layer " + std::to_string(j);
250 /* Use bins with a size of 10 cm per side. */
251 m_Spatial2DHitsEKLM[section - 1][j - 1] = new TH2F(name.c_str(), title.c_str(),
252 340 * 2 / 10, -340, 340,
253 340 * 2 / 10, -340, 340);
254 m_Spatial2DHitsEKLM[section - 1][j - 1]->GetXaxis()->SetTitle("X coordinate [cm]");
255 m_Spatial2DHitsEKLM[section - 1][j - 1]->GetYaxis()->SetTitle("Y coordinate [cm]");
256 }
257 }
258 }
259 oldDirectory->cd();
260}
261
263{
264 REG_HISTOGRAM;
265 m_RawFtsws.isOptional();
266 m_RawKlms.isOptional();
267 m_Digits.isOptional();
268 m_BklmHit1ds.isOptional();
269 m_Hit2ds.isOptional();
270}
271
273{
274 /* DAQ inclusion. */
275 m_DAQInclusion->Reset();
276 if (m_RawKlms.isValid())
277 m_DAQInclusion->Fill("Yes", 1);
278 else
279 m_DAQInclusion->Fill("No", 1);
280 /* Time. */
281 m_TimeRPC->Reset();
282 m_TimeScintillatorBKLM->Reset();
283 m_TimeScintillatorEKLM->Reset();
284 /* Plane hits. */
285 m_PlaneEKLM->Reset();
286 m_PlaneBKLMPhi->Reset();
287 m_PlaneBKLMZ->Reset();
288 /* Channel hits. */
290 for (KLMChannelIndex& klmSector : klmIndex) {
291 int nHistograms;
292 if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM)
293 nHistograms = m_ChannelHitHistogramsBKLM;
294 else
295 nHistograms = m_ChannelHitHistogramsEKLM;
296 KLMSectorNumber sector = klmSector.getKLMSectorNumber();
297 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sector);
298 for (int j = 0; j < nHistograms; j++)
299 m_ChannelHits[sectorIndex][j]->Reset();
300 }
301 /* Digits. */
302 m_DigitsKLM->Reset();
303 m_DigitsRPC->Reset();
306 m_DigitsMultiStripBKLM->Reset();
307 m_DigitsMultiStripEKLM->Reset();
308 /* Trigger bits. */
309 m_TriggerBitsBKLM->Reset();
310 m_TriggerBitsEKLM->Reset();
311 /* Injection information. */
312 m_DigitsAfterLERInj->Reset();
313 m_TriggersLERInj->Reset();
314 m_DigitsAfterHERInj->Reset();
315 m_TriggersHERInj->Reset();
316 /* Spatial 2D hits distributions. */
318 for (KLMChannelIndex& klmSection : klmIndex) {
319 KLMSubdetectorNumber subdetector = klmSection.getSubdetector();
320 if (subdetector == KLMElementNumbers::c_EKLM) {
321 KLMSectionNumber section = klmSection.getSection();
322 int maximalLayerNumber =
324 for (int j = 1; j <= maximalLayerNumber; ++j)
325 m_Spatial2DHitsEKLM[section - 1][j - 1]->Reset();
326 }
327 }
328}
329
331{
332 int nDigits = m_Digits.getEntries();
333 int nDigitsRPC = 0, nDigitsScintillatorBKLM = 0, nDigitsScintillatorEKLM = 0;
334 int nDigitsMultiStripBKLM = 0, nDigitsMultiStripEKLM = 0;
335 for (const KLMDigit& digit : m_Digits) {
336 /*
337 * Reject digits that are below the threshold (such digits may appear
338 * for simulated events).
339 */
340 if (!digit.isGood())
341 continue;
342 if (digit.getSubdetector() == KLMElementNumbers::c_EKLM) {
343 nDigitsScintillatorEKLM++;
344 int section = digit.getSection();
345 int sector = digit.getSector();
346 int layer = digit.getLayer();
347 int plane = digit.getPlane();
348 int strip = digit.getStrip();
349 if (not digit.isMultiStrip()) {
350 KLMSectorNumber klmSector = m_ElementNumbers->sectorNumberEKLM(section, sector);
351 KLMSectorNumber klmSectorIndex = m_SectorArrayIndex->getIndex(klmSector);
352 KLMChannelNumber channel = m_ElementNumbers->channelNumberEKLM(section, sector, layer, plane, strip);
353 KLMChannelNumber channelIndex = m_ChannelArrayIndex->getIndex(channel);
354 for (int j = 0; j < m_ChannelHitHistogramsEKLM; j++) {
355 double xMin = m_ChannelHits[klmSectorIndex][j]->GetXaxis()->GetXmin();
356 double xMax = m_ChannelHits[klmSectorIndex][j]->GetXaxis()->GetXmax();
357 if ((xMin > channelIndex) || (xMax < channelIndex))
358 continue;
359 m_ChannelHits[klmSectorIndex][j]->Fill(channelIndex);
360 }
361 } else
362 nDigitsMultiStripEKLM++;
363 int planeGlobal = m_eklmElementNumbers->planeNumber(section, layer, sector, plane);
364 m_PlaneEKLM->Fill(planeGlobal);
365 m_TimeScintillatorEKLM->Fill(digit.getTime());
366 if (digit.isMultiStrip()) {
367 KLMDigitRaw* digitRaw = digit.getRelated<KLMDigitRaw>();
368 if (digitRaw) {
369 uint16_t triggerBits = digitRaw->getTriggerBits();
370 if ((triggerBits & 0x1) != 0)
372 if ((triggerBits & 0x2) != 0)
374 if ((triggerBits & 0x4) != 0)
376 if ((triggerBits & 0x8) != 0)
378 }
379 }
380 } else if (digit.getSubdetector() == KLMElementNumbers::c_BKLM) {
381 int section = digit.getSection();
382 int sector = digit.getSector();
383 int layer = digit.getLayer();
384 int plane = digit.getPlane();
385 int strip = digit.getStrip();
386 if (not digit.isMultiStrip()) {
387 KLMSectorNumber klmSector = m_ElementNumbers->sectorNumberBKLM(section, sector);
388 KLMSectorNumber klmSectorIndex = m_SectorArrayIndex->getIndex(klmSector);
389 KLMChannelNumber channel = m_ElementNumbers->channelNumberBKLM(section, sector, layer, plane, strip);
390 KLMChannelNumber channelIndex = m_ChannelArrayIndex->getIndex(channel);
391 for (int j = 0; j < m_ChannelHitHistogramsBKLM; j++) {
392 double xMin = m_ChannelHits[klmSectorIndex][j]->GetXaxis()->GetXmin();
393 double xMax = m_ChannelHits[klmSectorIndex][j]->GetXaxis()->GetXmax();
394 if ((xMin > channelIndex) || (xMax < channelIndex))
395 continue;
396 m_ChannelHits[klmSectorIndex][j]->Fill(channelIndex);
397 }
398 } else
399 nDigitsMultiStripBKLM++;
400 if (digit.inRPC()) {
401 nDigitsRPC++;
402 m_TimeRPC->Fill(digit.getTime());
403 } else {
404 nDigitsScintillatorBKLM++;
405 m_TimeScintillatorBKLM->Fill(digit.getTime());
406 }
407 if (digit.isMultiStrip()) {
408 KLMDigitRaw* digitRaw = digit.getRelated<KLMDigitRaw>();
409 if (digitRaw) {
410 uint16_t triggerBits = digitRaw->getTriggerBits();
411 if ((triggerBits & 0x1) != 0)
413 if ((triggerBits & 0x2) != 0)
415 if ((triggerBits & 0x4) != 0)
417 if ((triggerBits & 0x8) != 0)
419 }
420 }
421 } else
422 B2FATAL("Not a BKLM or a EKLM digit, something went really wrong.");
423 }
424 for (const BKLMHit1d& hit1d : m_BklmHit1ds) {
425 int section = hit1d.getSection();
426 int sector = hit1d.getSector();
427 int layer = hit1d.getLayer();
429 section, sector, layer);
430 if (hit1d.isPhiReadout())
431 m_PlaneBKLMPhi->Fill(layerGlobal);
432 else
433 m_PlaneBKLMZ->Fill(layerGlobal);
434 }
435 /* Digits. */
436 m_DigitsKLM->Fill((double)nDigits);
437 m_DigitsRPC->Fill((double)nDigitsRPC);
438 m_DigitsScintillatorBKLM->Fill((double)nDigitsScintillatorBKLM);
439 m_DigitsScintillatorEKLM->Fill((double)nDigitsScintillatorEKLM);
440 if (nDigitsMultiStripBKLM > 0)
441 m_DigitsMultiStripBKLM->Fill((double)nDigitsMultiStripBKLM);
442 if (nDigitsMultiStripEKLM > 0)
443 m_DigitsMultiStripEKLM->Fill((double)nDigitsMultiStripEKLM);
444 /* Injection information. */
445 for (RawFTSW& rawFtsw : m_RawFtsws) {
446 unsigned int difference = rawFtsw.GetTimeSinceLastInjection(0);
447 if (difference != 0x7FFFFFFF) {
448 /* 127 MHz clock ticks to us, inexact rounding. */
449 float differenceInUs = difference / 127.;
450 if (rawFtsw.GetIsHER(0)) {
451 m_DigitsAfterHERInj->Fill(differenceInUs, nDigits);
452 m_TriggersHERInj->Fill(differenceInUs);
453 } else {
454 m_DigitsAfterLERInj->Fill(differenceInUs, nDigits);
455 m_TriggersLERInj->Fill(differenceInUs);
456 }
457 }
458 /*
459 * Usually, only one RawFTSW object is stored per event.
460 * If there are more, ignore the others.
461 */
462 break;
463 }
464 /* Spatial 2D hits distributions. */
465 for (const KLMHit2d& hit2d : m_Hit2ds) {
466 if (hit2d.getSubdetector() != KLMElementNumbers::c_EKLM)
467 continue;
468 int section = hit2d.getSection();
469 int layer = hit2d.getLayer();
470 m_Spatial2DHitsEKLM[section - 1][layer - 1]->Fill(hit2d.getPositionX(), hit2d.getPositionY());
471 }
472}
473
475{
476}
477
479{
480}
static constexpr int getMaximalSectorGlobalNumber()
Get maximal sector global number.
static int layerGlobalNumber(int section, int sector, int layer)
Get layer global number.
Store one reconstructed BKLM 1D hit as a ROOT object.
Definition: BKLMHit1d.h:30
EKLM element numbers.
int getMaximalDetectorLayerNumber(int section) const
Get maximal detector layer number.
int planeNumber(int section, int layer, int sector, int plane) const
Get plane number.
static constexpr int getMaximalSectorGlobalNumberKLMOrder()
Get maximal sector global number with KLM ordering (section, sector).
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
KLM channel array index.
KLM channel index.
int getSection() const
Get section.
int getSector() const
Get sector.
KLMChannelNumber getKLMChannelNumber() const
Get KLM channel number.
void setIndexLevel(enum IndexLevel indexLevel)
Set index level.
double m_BKLMScintTimeMin
Min time for BKLM Scint.
Definition: KLMDQMModule.h:116
TH1F * m_DigitsMultiStripEKLM
Number of multi-strip digits: EKLM scintillators.
Definition: KLMDQMModule.h:175
StoreArray< RawFTSW > m_RawFtsws
Raw FTSW.
Definition: KLMDQMModule.h:211
const int m_ChannelHitHistogramsEKLM
Number of channel hit histograms per sector for EKLM.
Definition: KLMDQMModule.h:107
TH1F * m_PlaneBKLMPhi
Plane occupancy: BKLM, phi readout.
Definition: KLMDQMModule.h:143
TH1F * m_TriggerBitsBKLM
Trigger bits: BKLM scintillators.
Definition: KLMDQMModule.h:178
~KLMDQMModule()
Destructor.
Definition: KLMDQMModule.cc:69
TH1F * m_TimeScintillatorEKLM
Time: EKLM scintillators.
Definition: KLMDQMModule.h:140
double m_RPCTimeMax
Max time for RPC.
Definition: KLMDQMModule.h:113
StoreArray< KLMDigit > m_Digits
KLM digits.
Definition: KLMDQMModule.h:217
double m_RPCTimeMin
Min time for RPC.
Definition: KLMDQMModule.h:110
void initialize() override
Initializer.
TH1F * m_DigitsScintillatorEKLM
Number of digits: EKLM scintillators.
Definition: KLMDQMModule.h:169
TH1F ** m_ChannelHits[EKLMElementNumbers::getMaximalSectorGlobalNumberKLMOrder()+BKLMElementNumbers::getMaximalSectorGlobalNumber()]
Number of hits per channel.
Definition: KLMDQMModule.h:154
TH1F * m_TimeRPC
Time: BKLM RPCs.
Definition: KLMDQMModule.h:134
TH1F * m_DigitsKLM
Number of digits: whole KLM.
Definition: KLMDQMModule.h:160
void event() override
This method is called for each event.
const KLMElementNumbers * m_ElementNumbers
KLM element numbers.
Definition: KLMDQMModule.h:205
KLMDQMModule()
Constructor.
Definition: KLMDQMModule.cc:23
double m_EKLMScintTimeMax
Max time for EKLM Scint.
Definition: KLMDQMModule.h:125
const int m_ChannelHitHistogramsBKLM
Number of channel hit histograms per sector for BKLM.
Definition: KLMDQMModule.h:104
TH1F * m_TriggerBitsEKLM
Trigger bits: EKLM scintillators.
Definition: KLMDQMModule.h:181
void endRun() override
This method is called if the current run ends.
const EKLMElementNumbers * m_eklmElementNumbers
Element numbers.
Definition: KLMDQMModule.h:208
double m_BKLMScintTimeMax
Max time for BKLM Scint.
Definition: KLMDQMModule.h:119
void terminate() override
This method is called at the end of the event processing.
TH1F * m_PlaneBKLMZ
Plane occupancy: BKLM, z readout.
Definition: KLMDQMModule.h:146
TH1F * m_TimeScintillatorBKLM
Time: BKLM scintillators.
Definition: KLMDQMModule.h:137
TH1F * m_DigitsRPC
Number of digits: BKLM RPCs.
Definition: KLMDQMModule.h:163
TH1F * m_TriggersHERInj
Histogram to be used for normalization of occupancy after HER injection.
Definition: KLMDQMModule.h:193
void beginRun() override
Called when entering a new run.
TH1F * m_DigitsAfterHERInj
Number of KLM Digits after LER injection.
Definition: KLMDQMModule.h:190
TH1F * m_DigitsScintillatorBKLM
Number of digits: BKLM scintillators.
Definition: KLMDQMModule.h:166
TH1F * m_TriggersLERInj
Histogram to be used for normalization of occupancy after LER injection.
Definition: KLMDQMModule.h:187
TH1F * m_DigitsMultiStripBKLM
Number of multi-strip digits: BKLM scintillators.
Definition: KLMDQMModule.h:172
const KLMSectorArrayIndex * m_SectorArrayIndex
KLM sector array index.
Definition: KLMDQMModule.h:202
StoreArray< BKLMHit1d > m_BklmHit1ds
BKLM 1d hits.
Definition: KLMDQMModule.h:220
TH1F * m_DigitsAfterLERInj
Number of KLM Digits after LER injection.
Definition: KLMDQMModule.h:184
TH1F * m_DAQInclusion
KLM DAQ inclusion.
Definition: KLMDQMModule.h:131
StoreArray< RawKLM > m_RawKlms
Raw KLM.
Definition: KLMDQMModule.h:214
TH1F * m_PlaneEKLM
Plane occupancy: EKLM.
Definition: KLMDQMModule.h:149
TH2F ** m_Spatial2DHitsEKLM[EKLMElementNumbers::getMaximalSectionNumber()]
Spatial distribution of EKLM 2d hits per layer.
Definition: KLMDQMModule.h:196
TH1F * m_MaskedChannelsPerSector
Masked channels per sector.
Definition: KLMDQMModule.h:157
double m_EKLMScintTimeMin
Min time for EKLM Scint.
Definition: KLMDQMModule.h:122
const KLMChannelArrayIndex * m_ChannelArrayIndex
KLM channel array index.
Definition: KLMDQMModule.h:199
std::string m_HistogramDirectoryName
Directory for KLM DQM histograms in ROOT file.
Definition: KLMDQMModule.h:128
void defineHisto() override
Definition of the histograms.
Definition: KLMDQMModule.cc:89
StoreArray< KLMHit2d > m_Hit2ds
KLM 2d hits.
Definition: KLMDQMModule.h:223
Class to store the raw words from the unpacker, digit-by-digit.
Definition: KLMDigitRaw.h:29
uint16_t getTriggerBits()
Get trigger bits.
Definition: KLMDigitRaw.h:159
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:29
uint16_t getNElements() const
Get number of elements.
uint16_t getIndex(uint16_t number) const
Get element index.
KLM element numbers.
KLMSectorNumber sectorNumberEKLM(int section, int sector) const
Get sector number for EKLM.
KLMChannelNumber channelNumberBKLM(int section, int sector, int layer, int plane, int strip) const
Get channel number for BKLM.
KLMSectorNumber sectorNumberBKLM(int section, int sector) const
Get sector number for BKLM.
std::string getSectorDAQName(int subdetector, int section, int sector) const
Get DAQ name for a given sector.
KLMChannelNumber channelNumberEKLM(int section, int sector, int layer, int plane, int strip) const
Get channel number for EKLM.
KLM 2d hit.
Definition: KLMHit2d.h:33
KLM sector array index.
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
The Raw FTSW class.
Definition: RawFTSW.h:30
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
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
uint16_t KLMSectorNumber
Sector number.
uint16_t KLMChannelNumber
Channel number.
uint16_t KLMSubdetectorNumber
Subdetector number.
uint16_t KLMSectionNumber
Section number.
Abstract base class for different kinds of events.