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_FE_BKLM_Layer_0(nullptr),
46 m_FE_BKLM_Layer_1(nullptr),
47 m_FE_EKLM_Plane_0(nullptr),
48 m_FE_EKLM_Plane_1(nullptr),
49 m_ChannelArrayIndex{&(KLMChannelArrayIndex::Instance())},
50 m_SectorArrayIndex{&(KLMSectorArrayIndex::Instance())},
51 m_ElementNumbers{&(KLMElementNumbers::Instance())},
52 m_eklmElementNumbers{&(EKLMElementNumbers::Instance())}
53{
54 setDescription("KLM data quality monitor.");
56 addParam("histogramDirectoryName", m_HistogramDirectoryName,
57 "Directory for KLM DQM histograms in ROOT file.",
58 std::string("KLM"));
59 addParam("RPCTimeMin", m_RPCTimeMin,
60 "Min time for RPC time histogram.", double(-1223.5));
61 addParam("RPCTimeMax", m_RPCTimeMax,
62 "Max time for RPC time histogram.", double(-199.5));
63 addParam("BKLMScintTimeMin", m_BKLMScintTimeMin,
64 "Min time for BKLM scintillator time histogram.", double(-5300));
65 addParam("BKLMScintTimeMax", m_BKLMScintTimeMax,
66 "Max time for BKLM scintillator time histogram.", double(-4300));
67 addParam("EKLMScintTimeMin", m_EKLMScintTimeMin,
68 "Min time for EKLM scintillator time histogram.", double(-5300));
69 addParam("EKLMScintTimeMax", m_EKLMScintTimeMax,
70 "Max time for EKLM scintillator time histogram.", double(-4300));
71}
72
74{
76 for (KLMChannelIndex& klmSector : klmIndex) {
77 KLMSectorNumber sector = klmSector.getKLMSectorNumber();
78 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sector);
79 if (m_ChannelHits[sectorIndex] != nullptr)
80 delete[] m_ChannelHits[sectorIndex];
81 }
83 for (KLMChannelIndex& klmSection : klmIndex) {
84 KLMSubdetectorNumber subdetector = klmSection.getSubdetector();
85 if (subdetector == KLMElementNumbers::c_EKLM) {
86 KLMSectionNumber section = klmSection.getSection();
87 if (m_Spatial2DHitsEKLM[section - 1] != nullptr)
88 delete[] m_Spatial2DHitsEKLM[section - 1];
89 }
90 }
91}
92
94{
95 TDirectory* oldDirectory, *newDirectory;
96 oldDirectory = gDirectory;
97 newDirectory = oldDirectory->mkdir(m_HistogramDirectoryName.c_str());
98 newDirectory->cd();
99 /* DAQ inclusion. */
100 m_DAQInclusion = new TH1F("daq_inclusion", "Is KLM included in DAQ?", 2, 0.0, 2.0);
101 m_DAQInclusion->GetXaxis()->SetBinLabel(1, "No");
102 m_DAQInclusion->GetXaxis()->SetBinLabel(2, "Yes");
103 /* Time histograms. */
104 m_TimeRPC = new TH1F("time_rpc", "RPC hit time", 128, m_RPCTimeMin, m_RPCTimeMax);
105 m_TimeRPC->GetXaxis()->SetTitle("Time, ns");
107 new TH1F("time_scintillator_bklm", "Scintillator hit time (BKLM)",
109 m_TimeScintillatorBKLM->GetXaxis()->SetTitle("Time, ns");
111 new TH1F("time_scintillator_eklm", "Scintillator hit time (EKLM)",
113 m_TimeScintillatorEKLM->GetXaxis()->SetTitle("Time, ns");
114 /* Number of hits per plane. */
115 m_PlaneBKLMPhi = new TH1F("plane_bklm_phi",
116 "BKLM plane occupancy (#phi readout)",
117 240, 0.5, 240.5);
118 m_PlaneBKLMPhi->GetXaxis()->SetTitle("Layer number");
119 m_PlaneBKLMZ = new TH1F("plane_bklm_z",
120 "BKLM plane occupancy (z readout)",
121 240, 0.5, 240.5);
122 m_PlaneBKLMZ->GetXaxis()->SetTitle("Layer number");
123 m_PlaneEKLM = new TH1F("plane_eklm", "EKLM plane occupancy (both readouts)", 208, 0.5, 208.5);
124 m_PlaneEKLM->GetXaxis()->SetTitle("Plane number");
125 /* Number of hits per channel. */
126 int nChannelHistograms =
131 KLMChannelNumber* firstChannelNumbers =
132 new KLMChannelNumber[nChannelHistograms + 1];
133 int i = 0;
135 for (KLMChannelIndex& klmSector : klmIndex) {
136 KLMChannelIndex klmChannel(klmSector);
138 KLMChannelNumber channel = klmChannel.getKLMChannelNumber();
139 firstChannelNumbers[i] = m_ChannelArrayIndex->getIndex(channel);
140 if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM) {
142 klmChannel.getSection(), klmChannel.getSector(), 8, 0, 1);
143 firstChannelNumbers[i + 1] = m_ChannelArrayIndex->getIndex(channel);
144 i += 2;
145 } else {
146 int layerIncrease = (klmSector.getSection() == 1) ? 4 : 5;
148 klmChannel.getSection(), klmChannel.getSector(),
149 1 + layerIncrease, 1, 1);
150 firstChannelNumbers[i + 1] = m_ChannelArrayIndex->getIndex(channel);
152 klmChannel.getSection(), klmChannel.getSector(),
153 1 + layerIncrease * 2, 1, 1);
154 firstChannelNumbers[i + 2] = m_ChannelArrayIndex->getIndex(channel);
155 i += 3;
156 }
157 }
158 firstChannelNumbers[nChannelHistograms] = m_ChannelArrayIndex->getNElements();
159 i = 0;
161 for (KLMChannelIndex& klmSector : klmIndex) {
162 int nHistograms;
163 if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM)
164 nHistograms = m_ChannelHitHistogramsBKLM;
165 else
166 nHistograms = m_ChannelHitHistogramsEKLM;
167 KLMSectorNumber sector = klmSector.getKLMSectorNumber();
168 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sector);
169 m_ChannelHits[sectorIndex] = new TH1F*[nHistograms];
170 for (int j = 0; j < nHistograms; j++) {
171 std::string name =
172 "strip_hits_subdetector_" + std::to_string(klmSector.getSubdetector()) +
173 "_section_" + std::to_string(klmSector.getSection()) +
174 "_sector_" + std::to_string(klmSector.getSector()) +
175 "_" + std::to_string(j);
176 std::string title = "Sector " + std::to_string(klmSector.getSector()) + " -- " +
177 m_ElementNumbers->getSectorDAQName(klmSector.getSubdetector(), klmSector.getSection(), klmSector.getSector());
178 m_ChannelHits[sectorIndex][j] = new TH1F(
179 name.c_str(), title.c_str(),
180 firstChannelNumbers[i + 1] - firstChannelNumbers[i],
181 firstChannelNumbers[i] - 0.5, firstChannelNumbers[i + 1] - 0.5);
182 i++;
183 }
184 }
185 delete[] firstChannelNumbers;
186 /* Masked channels per sector:
187 * it is defined here, but filled by the analysis module. */
189 m_MaskedChannelsPerSector = new TH1F("masked_channels", "Number of masked channels per sector",
190 totalSectors, -0.5, totalSectors - 0.5);
192 for (KLMChannelIndex& klmSector : klmIndex) {
193 std::string label = m_ElementNumbers->getSectorDAQName(klmSector.getSubdetector(), klmSector.getSection(), klmSector.getSector());
194 KLMSectorNumber sector = klmSector.getKLMSectorNumber();
195 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sector);
196 m_MaskedChannelsPerSector->GetXaxis()->SetBinLabel(sectorIndex + 1, label.c_str());
197 }
198 /* Number of digits. */
199 m_DigitsKLM = new TH1F("digits_klm", "Number of KLM digits",
200 250.0, 0.0, 250.0);
201 m_DigitsKLM->GetXaxis()->SetTitle("Number of digits");
202 m_DigitsRPC = new TH1F("digits_rpc", "Number of RPC digits",
203 250.0, 0.0, 250.0);
204 m_DigitsRPC->GetXaxis()->SetTitle("Number of digits");
205 m_DigitsScintillatorBKLM = new TH1F("digits_scintillator_bklm", "Number of scintillator digits (BKLM)",
206 250.0, 0.0, 250.0);
207 m_DigitsScintillatorBKLM->GetXaxis()->SetTitle("Number of digits");
208 m_DigitsScintillatorEKLM = new TH1F("digits_scintillator_eklm", "Number of scintillator digits (EKLM)",
209 250.0, 0.0, 250.0);
210 m_DigitsScintillatorEKLM->GetXaxis()->SetTitle("Number of digits");
211 m_DigitsMultiStripBKLM = new TH1F("digits_multi_bklm", "Number of multi-strip digits (BKLM)",
212 50.0, 0.0, 50.0);
213 m_DigitsMultiStripBKLM->GetXaxis()->SetTitle("Number of multi-strip digits");
214 m_DigitsMultiStripEKLM = new TH1F("digits_multi_eklm", "Number of multi-strip digits (EKLM)",
215 50.0, 0.0, 50.0);
216 m_DigitsMultiStripEKLM->GetXaxis()->SetTitle("Number of multi-strip digits");
217 /* Trigger bits. */
218 m_TriggerBitsBKLM = new TH1F("trigger_bits_bklm", "Trigger bits of multi-strip digits (BKLM)",
219 (double)c_0x1, (double)c_0x8, (double)c_0x1 + 1.0);
220 m_TriggerBitsBKLM->GetXaxis()->SetBinLabel(c_0x8, "0x8");
221 m_TriggerBitsBKLM->GetXaxis()->SetBinLabel(c_0x4, "0x4");
222 m_TriggerBitsBKLM->GetXaxis()->SetBinLabel(c_0x2, "0x2");
223 m_TriggerBitsBKLM->GetXaxis()->SetBinLabel(c_0x1, "0x1");
224 m_TriggerBitsEKLM = new TH1F("trigger_bits_eklm", "Trigger bits of multi-strip digits (EKLM)",
225 (double)c_0x1, (double)c_0x8, (double)c_0x1 + 1.0);
226 m_TriggerBitsEKLM->GetXaxis()->SetBinLabel(c_0x8, "0x8");
227 m_TriggerBitsEKLM->GetXaxis()->SetBinLabel(c_0x4, "0x4");
228 m_TriggerBitsEKLM->GetXaxis()->SetBinLabel(c_0x2, "0x2");
229 m_TriggerBitsEKLM->GetXaxis()->SetBinLabel(c_0x1, "0x1");
230 /* Number of digits after injection */
231 /* For the histograms below, we use the same style as for other subdetectors. */
232 m_DigitsAfterLERInj = new TH1F("KLMOccInjLER", "KLM digits after LER injection / Time;Time [#mus];Number of KLM digits / (5 #mus)",
233 4000, 0, 20000);
234 m_TriggersLERInj = new TH1F("KLMTrigInjLER", "Triggers after KER injection / Time;Time [#mus];Number of triggers / (5 #mus)",
235 4000, 0, 20000);
236 m_DigitsAfterHERInj = new TH1F("KLMOccInjHER", "KLM digits after HER injection / Time;Time [#mus];Number of KLM digits / (5 #mus)",
237 4000, 0, 20000);
238 m_TriggersHERInj = new TH1F("KLMTrigInjHER", "Triggers after HER injection / Time;Time [#mus];Number of triggers / (5 #mus)",
239 4000, 0, 20000);
240 /* Spatial distribution of EKLM 2d hits per layer. */
242 for (KLMChannelIndex& klmSection : klmIndex) {
243 KLMSubdetectorNumber subdetector = klmSection.getSubdetector();
244 if (subdetector == KLMElementNumbers::c_EKLM) {
245 KLMSectionNumber section = klmSection.getSection();
246 int maximalLayerNumber = m_eklmElementNumbers->getMaximalDetectorLayerNumber(section);
247 m_Spatial2DHitsEKLM[section - 1] = new TH2F*[maximalLayerNumber];
248 std::string sectionName = (section == EKLMElementNumbers::c_ForwardSection) ? "Forward" : "Backward";
249 for (int j = 1; j <= maximalLayerNumber; ++j) {
250 std::string name = "spatial_2d_hits_subdetector_" + std::to_string(subdetector) +
251 "_section_" + std::to_string(section) +
252 "_layer_" + std::to_string(j);
253 std::string title = "Endcap " + sectionName + " , Layer " + std::to_string(j);
254 /* Use bins with a size of 10 cm per side. */
255 m_Spatial2DHitsEKLM[section - 1][j - 1] = new TH2F(name.c_str(), title.c_str(),
256 340 * 2 / 10, -340, 340,
257 340 * 2 / 10, -340, 340);
258 m_Spatial2DHitsEKLM[section - 1][j - 1]->GetXaxis()->SetTitle("X coordinate [cm]");
259 m_Spatial2DHitsEKLM[section - 1][j - 1]->GetYaxis()->SetTitle("Y coordinate [cm]");
260 }
261 }
262 }
263 // Feature extraction status histogram for BKLM
264 int bklmSectors = BKLMElementNumbers::getMaximalSectorGlobalNumber(); // 16
265 int eklmPlanes = EKLMElementNumbers::getMaximalPlaneGlobalNumber(); // 208
266
267 m_FE_BKLM_Layer_0 = new TH1F("feStatus_bklm_scintillator_layers_0",
268 "BKLM Scintillator Standard Readout;FEE Card", bklmSectors * 2, 0.5, 0.5 + bklmSectors * 2);
269 m_FE_BKLM_Layer_1 = new TH1F("feStatus_bklm_scintillator_layers_1",
270 "BKLM Scintillator Feature Extraction;FEE Card", bklmSectors * 2, 0.5, 0.5 + bklmSectors * 2);
271 m_FE_EKLM_Plane_0 = new TH1F("feStatus_eklm_plane_0",
272 "EKLM Standard Readout;Plane number", eklmPlanes, 0.5, 0.5 + eklmPlanes);
273 m_FE_EKLM_Plane_1 = new TH1F("feStatus_eklm_plane_1",
274 "EKLM Feature Extraction;Plane number", eklmPlanes, 0.5, 0.5 + eklmPlanes);
275 oldDirectory->cd();
276}
277
279{
280 REG_HISTOGRAM;
281 m_RawFtsws.isOptional();
282 m_RawKlms.isOptional();
283 m_Digits.isOptional();
284 m_BklmHit1ds.isOptional();
285 m_Hit2ds.isOptional();
286}
287
289{
290 /* DAQ inclusion. */
291 m_DAQInclusion->Reset();
292 if (m_RawKlms.isValid())
293 m_DAQInclusion->Fill("Yes", 1);
294 else
295 m_DAQInclusion->Fill("No", 1);
296 /* Time. */
297 m_TimeRPC->Reset();
298 m_TimeScintillatorBKLM->Reset();
299 m_TimeScintillatorEKLM->Reset();
300 /* Plane hits. */
301 m_PlaneEKLM->Reset();
302 m_PlaneBKLMPhi->Reset();
303 m_PlaneBKLMZ->Reset();
304 /* Channel hits. */
306 for (KLMChannelIndex& klmSector : klmIndex) {
307 int nHistograms;
308 if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM)
309 nHistograms = m_ChannelHitHistogramsBKLM;
310 else
311 nHistograms = m_ChannelHitHistogramsEKLM;
312 KLMSectorNumber sector = klmSector.getKLMSectorNumber();
313 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sector);
314 for (int j = 0; j < nHistograms; j++)
315 m_ChannelHits[sectorIndex][j]->Reset();
316 }
317 /* Digits. */
318 m_DigitsKLM->Reset();
319 m_DigitsRPC->Reset();
322 m_DigitsMultiStripBKLM->Reset();
323 m_DigitsMultiStripEKLM->Reset();
324 /* Trigger bits. */
325 m_TriggerBitsBKLM->Reset();
326 m_TriggerBitsEKLM->Reset();
327 /* Injection information. */
328 m_DigitsAfterLERInj->Reset();
329 m_TriggersLERInj->Reset();
330 m_DigitsAfterHERInj->Reset();
331 m_TriggersHERInj->Reset();
332 /* Spatial 2D hits distributions. */
334 for (KLMChannelIndex& klmSection : klmIndex) {
335 KLMSubdetectorNumber subdetector = klmSection.getSubdetector();
336 if (subdetector == KLMElementNumbers::c_EKLM) {
337 KLMSectionNumber section = klmSection.getSection();
338 int maximalLayerNumber =
340 for (int j = 1; j <= maximalLayerNumber; ++j)
341 m_Spatial2DHitsEKLM[section - 1][j - 1]->Reset();
342 }
343 }
344 /* Feature extraction. */
345 m_FE_BKLM_Layer_0->Reset();
346 m_FE_BKLM_Layer_1->Reset();
347 m_FE_EKLM_Plane_0->Reset();
348 m_FE_EKLM_Plane_1->Reset();
349}
350
352{
353 int nDigits = m_Digits.getEntries();
354 int nDigitsRPC = 0, nDigitsScintillatorBKLM = 0, nDigitsScintillatorEKLM = 0;
355 int nDigitsMultiStripBKLM = 0, nDigitsMultiStripEKLM = 0;
356 for (const KLMDigit& digit : m_Digits) {
357 /*
358 * Reject digits that are below the threshold (such digits may appear
359 * for simulated events).
360 */
361
362 if (!digit.isGood())
363 continue;
364 KLMDigitRaw* digitRaw = digit.getRelated<KLMDigitRaw>();
365 if (digit.getSubdetector() == KLMElementNumbers::c_EKLM) {
366 nDigitsScintillatorEKLM++;
367 int section = digit.getSection();
368 int sector = digit.getSector();
369 int layer = digit.getLayer();
370 int plane = digit.getPlane();
371 int strip = digit.getStrip();
372 if (not digit.isMultiStrip()) {
373 KLMSectorNumber klmSector = m_ElementNumbers->sectorNumberEKLM(section, sector);
374 KLMSectorNumber klmSectorIndex = m_SectorArrayIndex->getIndex(klmSector);
375 KLMChannelNumber channel = m_ElementNumbers->channelNumberEKLM(section, sector, layer, plane, strip);
376 KLMChannelNumber channelIndex = m_ChannelArrayIndex->getIndex(channel);
377 for (int j = 0; j < m_ChannelHitHistogramsEKLM; j++) {
378 double xMin = m_ChannelHits[klmSectorIndex][j]->GetXaxis()->GetXmin();
379 double xMax = m_ChannelHits[klmSectorIndex][j]->GetXaxis()->GetXmax();
380 if ((xMin > channelIndex) || (xMax < channelIndex))
381 continue;
382 m_ChannelHits[klmSectorIndex][j]->Fill(channelIndex);
383 }
384 } else
385 nDigitsMultiStripEKLM++;
386 int planeGlobal = m_eklmElementNumbers->planeNumber(section, layer, sector, plane);
387 m_PlaneEKLM->Fill(planeGlobal);
388 m_TimeScintillatorEKLM->Fill(digit.getTime());
389 if (digit.isMultiStrip()) {
390 if (digitRaw) {
391 uint16_t triggerBits = digitRaw->getTriggerBits();
392 if ((triggerBits & 0x1) != 0)
394 if ((triggerBits & 0x2) != 0)
396 if ((triggerBits & 0x4) != 0)
398 if ((triggerBits & 0x8) != 0)
400 }
401 }
402 if (digitRaw) {
403 // Extract m_FE from m_word4
404 uint16_t feStatus = digitRaw->getFEStatus(); // Extract the most significant bit
405 if (feStatus != 0) {
406 m_FE_EKLM_Plane_1->Fill(planeGlobal);
407 } else {
408 m_FE_EKLM_Plane_0->Fill(planeGlobal);
409 }
410 }
411 } else if (digit.getSubdetector() == KLMElementNumbers::c_BKLM) {
412 int section = digit.getSection();
413 int sector = digit.getSector();
414 int layer = digit.getLayer();
415 int plane = digit.getPlane();
416 int strip = digit.getStrip();
417
418 KLMSectorNumber klmSector = m_ElementNumbers->sectorNumberBKLM(section, sector);
419 KLMSectorNumber klmSectorIndex = m_SectorArrayIndex->getIndex(klmSector);
420
421 if (not digit.isMultiStrip()) {
422 KLMChannelNumber channel = m_ElementNumbers->channelNumberBKLM(section, sector, layer, plane, strip);
423 KLMChannelNumber channelIndex = m_ChannelArrayIndex->getIndex(channel);
424 for (int j = 0; j < m_ChannelHitHistogramsBKLM; j++) {
425 double xMin = m_ChannelHits[klmSectorIndex][j]->GetXaxis()->GetXmin();
426 double xMax = m_ChannelHits[klmSectorIndex][j]->GetXaxis()->GetXmax();
427 if ((xMin > channelIndex) || (xMax < channelIndex))
428 continue;
429 m_ChannelHits[klmSectorIndex][j]->Fill(channelIndex);
430 }
431 } else
432 nDigitsMultiStripBKLM++;
433 if (digit.inRPC()) {
434 nDigitsRPC++;
435 m_TimeRPC->Fill(digit.getTime());
436 } else {
437 nDigitsScintillatorBKLM++;
438 m_TimeScintillatorBKLM->Fill(digit.getTime());
439 if (digitRaw) {
440 uint16_t feStatus = digitRaw->getFEStatus(); // Extract the most significant bit
441 if (feStatus != 0) {
442 m_FE_BKLM_Layer_1->Fill((klmSectorIndex) * 2 + layer);
443 } else {
444 m_FE_BKLM_Layer_0->Fill((klmSectorIndex) * 2 + layer);
445 }
446 }
447 }
448 if (digit.isMultiStrip()) {
449 if (digitRaw) {
450 uint16_t triggerBits = digitRaw->getTriggerBits();
451 if ((triggerBits & 0x1) != 0)
453 if ((triggerBits & 0x2) != 0)
455 if ((triggerBits & 0x4) != 0)
457 if ((triggerBits & 0x8) != 0)
459 }
460 }
461 } else
462 B2FATAL("Not a BKLM or a EKLM digit, something went really wrong.");
463 }
464 for (const BKLMHit1d& hit1d : m_BklmHit1ds) {
465 int section = hit1d.getSection();
466 int sector = hit1d.getSector();
467 int layer = hit1d.getLayer();
469 section, sector, layer);
470 if (hit1d.isPhiReadout())
471 m_PlaneBKLMPhi->Fill(layerGlobal);
472 else
473 m_PlaneBKLMZ->Fill(layerGlobal);
474 }
475 /* Digits. */
476 m_DigitsKLM->Fill((double)nDigits);
477 m_DigitsRPC->Fill((double)nDigitsRPC);
478 m_DigitsScintillatorBKLM->Fill((double)nDigitsScintillatorBKLM);
479 m_DigitsScintillatorEKLM->Fill((double)nDigitsScintillatorEKLM);
480 if (nDigitsMultiStripBKLM > 0)
481 m_DigitsMultiStripBKLM->Fill((double)nDigitsMultiStripBKLM);
482 if (nDigitsMultiStripEKLM > 0)
483 m_DigitsMultiStripEKLM->Fill((double)nDigitsMultiStripEKLM);
484 /* Injection information. */
485 for (RawFTSW& rawFtsw : m_RawFtsws) {
486 unsigned int difference = rawFtsw.GetTimeSinceLastInjection(0);
487 if (difference != 0x7FFFFFFF) {
488 /* 127 MHz clock ticks to us, inexact rounding. */
489 float differenceInUs = difference / 127.;
490 if (rawFtsw.GetIsHER(0)) {
491 m_DigitsAfterHERInj->Fill(differenceInUs, nDigits);
492 m_TriggersHERInj->Fill(differenceInUs);
493 } else {
494 m_DigitsAfterLERInj->Fill(differenceInUs, nDigits);
495 m_TriggersLERInj->Fill(differenceInUs);
496 }
497 }
498 /*
499 * Usually, only one RawFTSW object is stored per event.
500 * If there are more, ignore the others.
501 */
502 break;
503 }
504 /* Spatial 2D hits distributions. */
505 for (const KLMHit2d& hit2d : m_Hit2ds) {
506 if (hit2d.getSubdetector() != KLMElementNumbers::c_EKLM)
507 continue;
508 int section = hit2d.getSection();
509 int layer = hit2d.getLayer();
510 m_Spatial2DHitsEKLM[section - 1][layer - 1]->Fill(hit2d.getPositionX(), hit2d.getPositionY());
511 }
512}
513
515{
516}
517
519{
520}
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.
static constexpr int getMaximalPlaneGlobalNumber()
Get maximal plane global 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:223
const int m_ChannelHitHistogramsEKLM
Number of channel hit histograms per sector for EKLM.
Definition: KLMDQMModule.h:107
TH1F * m_FE_EKLM_Plane_1
feature extraction status for EKLM
Definition: KLMDQMModule.h:208
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:73
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:229
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:217
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:220
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_FE_BKLM_Layer_1
feature extraction status for BKLM Scintillator
Definition: KLMDQMModule.h:202
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:214
StoreArray< BKLMHit1d > m_BklmHit1ds
BKLM 1d hits.
Definition: KLMDQMModule.h:232
TH1F * m_DigitsAfterLERInj
Number of KLM Digits after LER injection.
Definition: KLMDQMModule.h:184
TH1F * m_DAQInclusion
KLM DAQ inclusion.
Definition: KLMDQMModule.h:131
TH1F * m_FE_EKLM_Plane_0
Standard Readout status for EKLM.
Definition: KLMDQMModule.h:205
StoreArray< RawKLM > m_RawKlms
Raw KLM.
Definition: KLMDQMModule.h:226
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:211
std::string m_HistogramDirectoryName
Directory for KLM DQM histograms in ROOT file.
Definition: KLMDQMModule.h:128
TH1F * m_FE_BKLM_Layer_0
Standard Readout status for BKLM Scintillator.
Definition: KLMDQMModule.h:199
void defineHisto() override
Definition of the histograms.
Definition: KLMDQMModule.cc:93
StoreArray< KLMHit2d > m_Hit2ds
KLM 2d hits.
Definition: KLMDQMModule.h:235
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
uint16_t getFEStatus()
Get FE.
Definition: KLMDigitRaw.h:175
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.