Belle II Software release-09-00-14
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_EventBackgroundTriggerSummary{nullptr},
42 m_DigitsAfterLERInj{nullptr},
43 m_TriggersLERInj{nullptr},
44 m_DigitsAfterHERInj{nullptr},
45 m_TriggersHERInj{nullptr},
46 m_FE_BKLM_Layer_0(nullptr),
47 m_FE_BKLM_Layer_1(nullptr),
48 m_FE_EKLM_Plane_0(nullptr),
49 m_FE_EKLM_Plane_1(nullptr),
50 m_ChannelArrayIndex{&(KLMChannelArrayIndex::Instance())},
51 m_SectorArrayIndex{&(KLMSectorArrayIndex::Instance())},
52 m_ElementNumbers{&(KLMElementNumbers::Instance())},
53 m_eklmElementNumbers{&(EKLMElementNumbers::Instance())},
54 m_klmTime{&(KLMTime::Instance())}
55{
56 setDescription("KLM data quality monitor.");
58 addParam("histogramDirectoryName", m_HistogramDirectoryName,
59 "Directory for KLM DQM histograms in ROOT file.",
60 std::string("KLM"));
61 addParam("RPCTimeMin", m_RPCTimeMin,
62 "Min time for RPC time histogram.", double(-1223.5));
63 addParam("RPCTimeMax", m_RPCTimeMax,
64 "Max time for RPC time histogram.", double(-199.5));
65 addParam("BKLMScintTimeMin", m_BKLMScintTimeMin,
66 "Min time for BKLM scintillator time histogram.", double(-5300));
67 addParam("BKLMScintTimeMax", m_BKLMScintTimeMax,
68 "Max time for BKLM scintillator time histogram.", double(-4300));
69 addParam("EKLMScintTimeMin", m_EKLMScintTimeMin,
70 "Min time for EKLM scintillator time histogram.", double(-5300));
71 addParam("EKLMScintTimeMax", m_EKLMScintTimeMax,
72 "Max time for EKLM scintillator time histogram.", double(-4300));
73 addParam("Revo9DCArrivalTimeMin", m_Revo9DCArrivalTimeMin,
74 "Min time for RPC hit revo9DCArrivalTime histogram.", double(-10000));
75 addParam("Revo9DCArrivalTimeMax", m_Revo9DCArrivalTimeMax,
76 "Max time for RPC hit revo9DCArrivalTime histogram.", double(3000));
77}
78
80{
82 for (KLMChannelIndex& klmSector : klmIndex) {
83 KLMSectorNumber sector = klmSector.getKLMSectorNumber();
84 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sector);
85 if (m_ChannelHits[sectorIndex] != nullptr)
86 delete[] m_ChannelHits[sectorIndex];
87 }
89 for (KLMChannelIndex& klmSection : klmIndex) {
90 KLMSubdetectorNumber subdetector = klmSection.getSubdetector();
91 if (subdetector == KLMElementNumbers::c_EKLM) {
92 KLMSectionNumber section = klmSection.getSection();
93 if (m_Spatial2DHitsEKLM[section - 1] != nullptr)
94 delete[] m_Spatial2DHitsEKLM[section - 1];
95 }
96 }
97}
98
100{
101 TDirectory* oldDirectory, *newDirectory;
102 oldDirectory = gDirectory;
103 newDirectory = oldDirectory->mkdir(m_HistogramDirectoryName.c_str());
104 newDirectory->cd();
105 /* DAQ inclusion. */
106 m_DAQInclusion = new TH1F("daq_inclusion", "Is KLM included in DAQ?", 2, 0.0, 2.0);
107 m_DAQInclusion->GetXaxis()->SetBinLabel(1, "No");
108 m_DAQInclusion->GetXaxis()->SetBinLabel(2, "Yes");
109 /* Time histograms. */
110 m_TimeRPC = new TH1F("time_rpc", "RPC hit time", 128, m_RPCTimeMin, m_RPCTimeMax);
111 m_TimeRPC->GetXaxis()->SetTitle("Time, ns");
113 new TH1F("time_scintillator_bklm", "Scintillator hit time (BKLM)",
115 m_TimeScintillatorBKLM->GetXaxis()->SetTitle("Time, ns");
117 new TH1F("time_scintillator_eklm", "Scintillator hit time (EKLM)",
119 m_TimeScintillatorEKLM->GetXaxis()->SetTitle("Time, ns");
120 m_TimeRevo9DCArrivalTime = new TH1F("time_revo9dc_arrival_time", "DC arrival hit time (RPC)",
122 m_TimeRevo9DCArrivalTime->GetXaxis()->SetTitle("Time, ns");
123 /* Number of hits per plane. */
124 m_PlaneBKLMPhi = new TH1F("plane_bklm_phi",
125 "BKLM plane occupancy (#phi readout)",
126 240, 0.5, 240.5);
127 m_PlaneBKLMPhi->GetXaxis()->SetTitle("Layer number");
128 m_PlaneBKLMZ = new TH1F("plane_bklm_z",
129 "BKLM plane occupancy (z readout)",
130 240, 0.5, 240.5);
131 m_PlaneBKLMZ->GetXaxis()->SetTitle("Layer number");
132 m_PlaneEKLM = new TH1F("plane_eklm", "EKLM plane occupancy (both readouts)", 208, 0.5, 208.5);
133 m_PlaneEKLM->GetXaxis()->SetTitle("Plane number");
134 /* Number of hits per channel. */
135 int nChannelHistograms =
140 KLMChannelNumber* firstChannelNumbers =
141 new KLMChannelNumber[nChannelHistograms + 1];
142 int i = 0;
144 for (KLMChannelIndex& klmSector : klmIndex) {
145 KLMChannelIndex klmChannel(klmSector);
147 KLMChannelNumber channel = klmChannel.getKLMChannelNumber();
148 firstChannelNumbers[i] = m_ChannelArrayIndex->getIndex(channel);
149 if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM) {
151 klmChannel.getSection(), klmChannel.getSector(), 8, 0, 1);
152 firstChannelNumbers[i + 1] = m_ChannelArrayIndex->getIndex(channel);
153 i += 2;
154 } else {
155 int layerIncrease = (klmSector.getSection() == 1) ? 4 : 5;
157 klmChannel.getSection(), klmChannel.getSector(),
158 1 + layerIncrease, 1, 1);
159 firstChannelNumbers[i + 1] = m_ChannelArrayIndex->getIndex(channel);
161 klmChannel.getSection(), klmChannel.getSector(),
162 1 + layerIncrease * 2, 1, 1);
163 firstChannelNumbers[i + 2] = m_ChannelArrayIndex->getIndex(channel);
164 i += 3;
165 }
166 }
167 firstChannelNumbers[nChannelHistograms] = m_ChannelArrayIndex->getNElements();
168 i = 0;
170 for (KLMChannelIndex& klmSector : klmIndex) {
171 int nHistograms;
172 if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM)
173 nHistograms = m_ChannelHitHistogramsBKLM;
174 else
175 nHistograms = m_ChannelHitHistogramsEKLM;
176 KLMSectorNumber sector = klmSector.getKLMSectorNumber();
177 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sector);
178 m_ChannelHits[sectorIndex] = new TH1F*[nHistograms];
179 for (int j = 0; j < nHistograms; j++) {
180 std::string name =
181 "strip_hits_subdetector_" + std::to_string(klmSector.getSubdetector()) +
182 "_section_" + std::to_string(klmSector.getSection()) +
183 "_sector_" + std::to_string(klmSector.getSector()) +
184 "_" + std::to_string(j);
185 std::string title = "Sector " + std::to_string(klmSector.getSector()) + " -- " +
186 m_ElementNumbers->getSectorDAQName(klmSector.getSubdetector(), klmSector.getSection(), klmSector.getSector());
187 m_ChannelHits[sectorIndex][j] = new TH1F(
188 name.c_str(), title.c_str(),
189 firstChannelNumbers[i + 1] - firstChannelNumbers[i],
190 firstChannelNumbers[i] - 0.5, firstChannelNumbers[i + 1] - 0.5);
191 i++;
192 }
193 }
194 delete[] firstChannelNumbers;
195 /* Masked channels per sector:
196 * it is defined here, but filled by the analysis module. */
198 m_MaskedChannelsPerSector = new TH1F("masked_channels", "Number of masked channels per sector",
199 totalSectors, -0.5, totalSectors - 0.5);
201 for (KLMChannelIndex& klmSector : klmIndex) {
202 std::string label = m_ElementNumbers->getSectorDAQName(klmSector.getSubdetector(), klmSector.getSection(), klmSector.getSector());
203 KLMSectorNumber sector = klmSector.getKLMSectorNumber();
204 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sector);
205 m_MaskedChannelsPerSector->GetXaxis()->SetBinLabel(sectorIndex + 1, label.c_str());
206 }
207 /* Number of digits. */
208 m_DigitsKLM = new TH1F("digits_klm", "Number of KLM digits",
209 250.0, 0.0, 250.0);
210 m_DigitsKLM->GetXaxis()->SetTitle("Number of digits");
211 m_DigitsRPC = new TH1F("digits_rpc", "Number of RPC digits",
212 250.0, 0.0, 250.0);
213 m_DigitsRPC->GetXaxis()->SetTitle("Number of digits");
214 m_DigitsScintillatorBKLM = new TH1F("digits_scintillator_bklm", "Number of scintillator digits (BKLM)",
215 250.0, 0.0, 250.0);
216 m_DigitsScintillatorBKLM->GetXaxis()->SetTitle("Number of digits");
217 m_DigitsScintillatorEKLM = new TH1F("digits_scintillator_eklm", "Number of scintillator digits (EKLM)",
218 250.0, 0.0, 250.0);
219 m_DigitsScintillatorEKLM->GetXaxis()->SetTitle("Number of digits");
220 m_DigitsMultiStripBKLM = new TH1F("digits_multi_bklm", "Number of multi-strip digits (BKLM)",
221 50.0, 0.0, 50.0);
222 m_DigitsMultiStripBKLM->GetXaxis()->SetTitle("Number of multi-strip digits");
223 m_DigitsMultiStripEKLM = new TH1F("digits_multi_eklm", "Number of multi-strip digits (EKLM)",
224 50.0, 0.0, 50.0);
225 m_DigitsMultiStripEKLM->GetXaxis()->SetTitle("Number of multi-strip digits");
226 /* Trigger bits. */
227 m_TriggerBitsBKLM = new TH1F("trigger_bits_bklm", "Trigger bits of multi-strip digits (BKLM)",
228 (double)c_0x1, (double)c_0x8, (double)c_0x1 + 1.0);
229 m_TriggerBitsBKLM->GetXaxis()->SetBinLabel(c_0x8, "0x8");
230 m_TriggerBitsBKLM->GetXaxis()->SetBinLabel(c_0x4, "0x4");
231 m_TriggerBitsBKLM->GetXaxis()->SetBinLabel(c_0x2, "0x2");
232 m_TriggerBitsBKLM->GetXaxis()->SetBinLabel(c_0x1, "0x1");
233 m_TriggerBitsEKLM = new TH1F("trigger_bits_eklm", "Trigger bits of multi-strip digits (EKLM)",
234 (double)c_0x1, (double)c_0x8, (double)c_0x1 + 1.0);
235 m_TriggerBitsEKLM->GetXaxis()->SetBinLabel(c_0x8, "0x8");
236 m_TriggerBitsEKLM->GetXaxis()->SetBinLabel(c_0x4, "0x4");
237 m_TriggerBitsEKLM->GetXaxis()->SetBinLabel(c_0x2, "0x2");
238 m_TriggerBitsEKLM->GetXaxis()->SetBinLabel(c_0x1, "0x1");
239
240 /* Event-level L1 trigger bits (TRGSummary). */
241 int nTimingBits = (int)c_KlmL1Triggers.size(); // TTYP_DPHY, TTYP_RAND, TTYP_POIS
242 m_EventBackgroundTriggerSummary = new TH1F("event_background_trigger_summary",
243 "Event background trigger summary;Trigger Decision;Events",
244 nTimingBits, 0.5, 0.5 + nTimingBits);
245 /* Keep labels explicit to match the actual trigger bits. */
246 m_EventBackgroundTriggerSummary->GetXaxis()->SetBinLabel(1, "TTYP_DPHY");
247 m_EventBackgroundTriggerSummary->GetXaxis()->SetBinLabel(2, "TTYP_RAND");
248 m_EventBackgroundTriggerSummary->GetXaxis()->SetBinLabel(3, "TTYP_POIS");
249
250 /* Number of digits after injection */
251 /* For the histograms below, we use the same style as for other subdetectors. */
252 m_DigitsAfterLERInj = new TH1F("KLMOccInjLER", "KLM digits after LER injection / Time;Time [#mus];Number of KLM digits / (5 #mus)",
253 4000, 0, 20000);
254 m_TriggersLERInj = new TH1F("KLMTrigInjLER", "Triggers after KER injection / Time;Time [#mus];Number of triggers / (5 #mus)",
255 4000, 0, 20000);
256 m_DigitsAfterHERInj = new TH1F("KLMOccInjHER", "KLM digits after HER injection / Time;Time [#mus];Number of KLM digits / (5 #mus)",
257 4000, 0, 20000);
258 m_TriggersHERInj = new TH1F("KLMTrigInjHER", "Triggers after HER injection / Time;Time [#mus];Number of triggers / (5 #mus)",
259 4000, 0, 20000);
260 /* Spatial distribution of EKLM 2d hits per layer. */
262 for (KLMChannelIndex& klmSection : klmIndex) {
263 KLMSubdetectorNumber subdetector = klmSection.getSubdetector();
264 if (subdetector == KLMElementNumbers::c_EKLM) {
265 KLMSectionNumber section = klmSection.getSection();
266 int maximalLayerNumber = m_eklmElementNumbers->getMaximalDetectorLayerNumber(section);
267 m_Spatial2DHitsEKLM[section - 1] = new TH2F*[maximalLayerNumber];
268 std::string sectionName = (section == EKLMElementNumbers::c_ForwardSection) ? "Forward" : "Backward";
269 for (int j = 1; j <= maximalLayerNumber; ++j) {
270 std::string name = "spatial_2d_hits_subdetector_" + std::to_string(subdetector) +
271 "_section_" + std::to_string(section) +
272 "_layer_" + std::to_string(j);
273 std::string title = "Endcap " + sectionName + " , Layer " + std::to_string(j);
274 /* Use bins with a size of 10 cm per side. */
275 m_Spatial2DHitsEKLM[section - 1][j - 1] = new TH2F(name.c_str(), title.c_str(),
276 340 * 2 / 10, -340, 340,
277 340 * 2 / 10, -340, 340);
278 m_Spatial2DHitsEKLM[section - 1][j - 1]->GetXaxis()->SetTitle("X coordinate [cm]");
279 m_Spatial2DHitsEKLM[section - 1][j - 1]->GetYaxis()->SetTitle("Y coordinate [cm]");
280 }
281 }
282 }
283 // Feature extraction status histogram for BKLM
284 int bklmSectors = BKLMElementNumbers::getMaximalSectorGlobalNumber(); // 16
285 int eklmPlanes = EKLMElementNumbers::getMaximalPlaneGlobalNumber(); // 208
286
287 m_FE_BKLM_Layer_0 = new TH1F("feStatus_bklm_scintillator_layers_0",
288 "BKLM Scintillator Standard Readout;FEE Card", bklmSectors * 2, 0.5, 0.5 + bklmSectors * 2);
289 m_FE_BKLM_Layer_1 = new TH1F("feStatus_bklm_scintillator_layers_1",
290 "BKLM Scintillator Feature Extraction;FEE Card", bklmSectors * 2, 0.5, 0.5 + bklmSectors * 2);
291 m_FE_EKLM_Plane_0 = new TH1F("feStatus_eklm_plane_0",
292 "EKLM Standard Readout;Plane number", eklmPlanes, 0.5, 0.5 + eklmPlanes);
293 m_FE_EKLM_Plane_1 = new TH1F("feStatus_eklm_plane_1",
294 "EKLM Feature Extraction;Plane number", eklmPlanes, 0.5, 0.5 + eklmPlanes);
295 oldDirectory->cd();
296}
297
299{
300 REG_HISTOGRAM;
301 m_RawFtsws.isOptional();
302 m_RawKlms.isOptional();
303 m_trgSummary.isOptional();
304 m_Digits.isOptional();
305 m_BklmHit1ds.isOptional();
306 m_Hit2ds.isOptional();
307}
308
310{
311 /* DAQ inclusion. */
312 m_DAQInclusion->Reset();
313 if (m_RawKlms.isValid())
314 m_DAQInclusion->Fill("Yes", 1);
315 else
316 m_DAQInclusion->Fill("No", 1);
317 /* Time. */
318 m_TimeRPC->Reset();
319 m_TimeScintillatorBKLM->Reset();
320 m_TimeScintillatorEKLM->Reset();
322 m_klmTime->updateConstants(); //to get correct CTime
323 /* Plane hits. */
324 m_PlaneEKLM->Reset();
325 m_PlaneBKLMPhi->Reset();
326 m_PlaneBKLMZ->Reset();
327 /* Channel hits. */
329 for (KLMChannelIndex& klmSector : klmIndex) {
330 int nHistograms;
331 if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM)
332 nHistograms = m_ChannelHitHistogramsBKLM;
333 else
334 nHistograms = m_ChannelHitHistogramsEKLM;
335 KLMSectorNumber sector = klmSector.getKLMSectorNumber();
336 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sector);
337 for (int j = 0; j < nHistograms; j++)
338 m_ChannelHits[sectorIndex][j]->Reset();
339 }
340 /* Digits. */
341 m_DigitsKLM->Reset();
342 m_DigitsRPC->Reset();
345 m_DigitsMultiStripBKLM->Reset();
346 m_DigitsMultiStripEKLM->Reset();
347 /* Trigger bits. */
348 m_TriggerBitsBKLM->Reset();
349 m_TriggerBitsEKLM->Reset();
350 /* Event-level background trigger summary. */
351 if (m_EventBackgroundTriggerSummary != nullptr)
353 /* Injection information. */
354 m_DigitsAfterLERInj->Reset();
355 m_TriggersLERInj->Reset();
356 m_DigitsAfterHERInj->Reset();
357 m_TriggersHERInj->Reset();
358 /* Spatial 2D hits distributions. */
360 for (KLMChannelIndex& klmSection : klmIndex) {
361 KLMSubdetectorNumber subdetector = klmSection.getSubdetector();
362 if (subdetector == KLMElementNumbers::c_EKLM) {
363 KLMSectionNumber section = klmSection.getSection();
364 int maximalLayerNumber =
366 for (int j = 1; j <= maximalLayerNumber; ++j)
367 m_Spatial2DHitsEKLM[section - 1][j - 1]->Reset();
368 }
369 }
370 /* Feature extraction. */
371 m_FE_BKLM_Layer_0->Reset();
372 m_FE_BKLM_Layer_1->Reset();
373 m_FE_EKLM_Plane_0->Reset();
374 m_FE_EKLM_Plane_1->Reset();
375}
376
378{
379 /* Event-level background trigger summary, filled once per event. */
380 if (m_trgSummary.isValid() && m_EventBackgroundTriggerSummary != nullptr) {
381 const int nTimingBits = (int)c_KlmL1Triggers.size();
382 // GDL Background triggers: TTYP_DPHY, TTYP_RAND, TTYP_POIS (event-level timing source)
383 for (int i = 0; i < nTimingBits; ++i) {
384 if (m_trgSummary->getTimType() == c_KlmL1Triggers[i])
385 m_EventBackgroundTriggerSummary->Fill((double)i + 1.0);
386 }
387 }
388
389 int nDigits = m_Digits.getEntries();
390 int nDigitsRPC = 0, nDigitsScintillatorBKLM = 0, nDigitsScintillatorEKLM = 0;
391 int nDigitsMultiStripBKLM = 0, nDigitsMultiStripEKLM = 0;
392 for (const KLMDigit& digit : m_Digits) {
393 /*
394 * Reject digits that are below the threshold (such digits may appear
395 * for simulated events).
396 */
397
398 if (!digit.isGood())
399 continue;
400 KLMDigitRaw* digitRaw = digit.getRelated<KLMDigitRaw>();
401 if (digit.getSubdetector() == KLMElementNumbers::c_EKLM) {
402 nDigitsScintillatorEKLM++;
403 int section = digit.getSection();
404 int sector = digit.getSector();
405 int layer = digit.getLayer();
406 int plane = digit.getPlane();
407 int strip = digit.getStrip();
408 if (not digit.isMultiStrip()) {
409 KLMSectorNumber klmSector = m_ElementNumbers->sectorNumberEKLM(section, sector);
410 KLMSectorNumber klmSectorIndex = m_SectorArrayIndex->getIndex(klmSector);
411 KLMChannelNumber channel = m_ElementNumbers->channelNumberEKLM(section, sector, layer, plane, strip);
412 KLMChannelNumber channelIndex = m_ChannelArrayIndex->getIndex(channel);
413 for (int j = 0; j < m_ChannelHitHistogramsEKLM; j++) {
414 double xMin = m_ChannelHits[klmSectorIndex][j]->GetXaxis()->GetXmin();
415 double xMax = m_ChannelHits[klmSectorIndex][j]->GetXaxis()->GetXmax();
416 if ((xMin > channelIndex) || (xMax < channelIndex))
417 continue;
418 m_ChannelHits[klmSectorIndex][j]->Fill(channelIndex);
419 }
420 } else
421 nDigitsMultiStripEKLM++;
422 int planeGlobal = m_eklmElementNumbers->planeNumber(section, layer, sector, plane);
423 m_PlaneEKLM->Fill(planeGlobal);
424 m_TimeScintillatorEKLM->Fill(digit.getTime());
425 if (digit.isMultiStrip()) {
426 if (digitRaw) {
427 uint16_t triggerBits = digitRaw->getTriggerBits();
428 if ((triggerBits & 0x1) != 0)
430 if ((triggerBits & 0x2) != 0)
432 if ((triggerBits & 0x4) != 0)
434 if ((triggerBits & 0x8) != 0)
436 }
437 }
438 if (digitRaw) {
439 // Extract m_FE from m_word4
440 uint16_t feStatus = digitRaw->getFEStatus(); // Extract the most significant bit
441 if (feStatus != 0) {
442 m_FE_EKLM_Plane_1->Fill(planeGlobal);
443 } else {
444 m_FE_EKLM_Plane_0->Fill(planeGlobal);
445 }
446 }
447 } else if (digit.getSubdetector() == KLMElementNumbers::c_BKLM) {
448 int section = digit.getSection();
449 int sector = digit.getSector();
450 int layer = digit.getLayer();
451 int plane = digit.getPlane();
452 int strip = digit.getStrip();
453
454 KLMSectorNumber klmSector = m_ElementNumbers->sectorNumberBKLM(section, sector);
455 KLMSectorNumber klmSectorIndex = m_SectorArrayIndex->getIndex(klmSector);
456
457 if (not digit.isMultiStrip()) {
458 KLMChannelNumber channel = m_ElementNumbers->channelNumberBKLM(section, sector, layer, plane, strip);
459 KLMChannelNumber channelIndex = m_ChannelArrayIndex->getIndex(channel);
460 for (int j = 0; j < m_ChannelHitHistogramsBKLM; j++) {
461 double xMin = m_ChannelHits[klmSectorIndex][j]->GetXaxis()->GetXmin();
462 double xMax = m_ChannelHits[klmSectorIndex][j]->GetXaxis()->GetXmax();
463 if ((xMin > channelIndex) || (xMax < channelIndex))
464 continue;
465 m_ChannelHits[klmSectorIndex][j]->Fill(channelIndex);
466 }
467 } else
468 nDigitsMultiStripBKLM++;
469 if (digit.inRPC()) {
470 nDigitsRPC++;
471 m_TimeRPC->Fill(digit.getTime());
472 m_TimeRevo9DCArrivalTime->Fill(digit.getRevo9DCArrivalTime() * m_klmTime->getCTimePeriod());
473 } else {
474 nDigitsScintillatorBKLM++;
475 m_TimeScintillatorBKLM->Fill(digit.getTime());
476 if (digitRaw) {
477 uint16_t feStatus = digitRaw->getFEStatus(); // Extract the most significant bit
478 if (feStatus != 0) {
479 m_FE_BKLM_Layer_1->Fill((klmSectorIndex) * 2 + layer);
480 } else {
481 m_FE_BKLM_Layer_0->Fill((klmSectorIndex) * 2 + layer);
482 }
483 }
484 }
485 if (digit.isMultiStrip()) {
486 if (digitRaw) {
487 uint16_t triggerBits = digitRaw->getTriggerBits();
488 if ((triggerBits & 0x1) != 0)
490 if ((triggerBits & 0x2) != 0)
492 if ((triggerBits & 0x4) != 0)
494 if ((triggerBits & 0x8) != 0)
496 }
497 }
498 } else
499 B2FATAL("Not a BKLM or a EKLM digit, something went really wrong.");
500 }
501 for (const BKLMHit1d& hit1d : m_BklmHit1ds) {
502 int section = hit1d.getSection();
503 int sector = hit1d.getSector();
504 int layer = hit1d.getLayer();
506 section, sector, layer);
507 if (hit1d.isPhiReadout())
508 m_PlaneBKLMPhi->Fill(layerGlobal);
509 else
510 m_PlaneBKLMZ->Fill(layerGlobal);
511 }
512 /* Digits. */
513 m_DigitsKLM->Fill((double)nDigits);
514 m_DigitsRPC->Fill((double)nDigitsRPC);
515 m_DigitsScintillatorBKLM->Fill((double)nDigitsScintillatorBKLM);
516 m_DigitsScintillatorEKLM->Fill((double)nDigitsScintillatorEKLM);
517 if (nDigitsMultiStripBKLM > 0)
518 m_DigitsMultiStripBKLM->Fill((double)nDigitsMultiStripBKLM);
519 if (nDigitsMultiStripEKLM > 0)
520 m_DigitsMultiStripEKLM->Fill((double)nDigitsMultiStripEKLM);
521 /* Injection information. */
522 for (RawFTSW& rawFtsw : m_RawFtsws) {
523 unsigned int difference = rawFtsw.GetTimeSinceLastInjection(0);
524 if (difference != 0x7FFFFFFF) {
525 /* 127 MHz clock ticks to us, inexact rounding. */
526 float differenceInUs = difference / 127.;
527 if (rawFtsw.GetIsHER(0)) {
528 m_DigitsAfterHERInj->Fill(differenceInUs, nDigits);
529 m_TriggersHERInj->Fill(differenceInUs);
530 } else {
531 m_DigitsAfterLERInj->Fill(differenceInUs, nDigits);
532 m_TriggersLERInj->Fill(differenceInUs);
533 }
534 }
535 /*
536 * Usually, only one RawFTSW object is stored per event.
537 * If there are more, ignore the others.
538 */
539 break;
540 }
541 /* Spatial 2D hits distributions. */
542 for (const KLMHit2d& hit2d : m_Hit2ds) {
543 if (hit2d.getSubdetector() != KLMElementNumbers::c_EKLM)
544 continue;
545 int section = hit2d.getSection();
546 int layer = hit2d.getLayer();
547 m_Spatial2DHitsEKLM[section - 1][layer - 1]->Fill(hit2d.getPositionX(), hit2d.getPositionY());
548 }
549}
550
552{
553}
554
556{
557}
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:124
TH1F * m_DigitsMultiStripEKLM
Number of multi-strip digits: EKLM scintillators.
Definition: KLMDQMModule.h:192
StoreArray< RawFTSW > m_RawFtsws
Raw FTSW.
Definition: KLMDQMModule.h:246
const int m_ChannelHitHistogramsEKLM
Number of channel hit histograms per sector for EKLM.
Definition: KLMDQMModule.h:115
TH1F * m_FE_EKLM_Plane_1
feature extraction status for EKLM
Definition: KLMDQMModule.h:228
TH1F * m_PlaneBKLMPhi
Plane occupancy: BKLM, phi readout.
Definition: KLMDQMModule.h:160
TH1F * m_TriggerBitsBKLM
Trigger bits: BKLM scintillators.
Definition: KLMDQMModule.h:195
~KLMDQMModule()
Destructor.
Definition: KLMDQMModule.cc:79
TH1F * m_TimeScintillatorEKLM
Time: EKLM scintillators.
Definition: KLMDQMModule.h:157
double m_RPCTimeMax
Max time for RPC.
Definition: KLMDQMModule.h:121
StoreArray< KLMDigit > m_Digits
KLM digits.
Definition: KLMDQMModule.h:255
double m_RPCTimeMin
Min time for RPC.
Definition: KLMDQMModule.h:118
void initialize() override
Initializer.
TH1F * m_DigitsScintillatorEKLM
Number of digits: EKLM scintillators.
Definition: KLMDQMModule.h:186
TH1F ** m_ChannelHits[EKLMElementNumbers::getMaximalSectorGlobalNumberKLMOrder()+BKLMElementNumbers::getMaximalSectorGlobalNumber()]
Number of hits per channel.
Definition: KLMDQMModule.h:171
TH1F * m_TimeRPC
Time: BKLM RPCs.
Definition: KLMDQMModule.h:148
TH1F * m_DigitsKLM
Number of digits: whole KLM.
Definition: KLMDQMModule.h:177
double m_Revo9DCArrivalTimeMax
Max time for revo9DCArrivalTime for RPC.
Definition: KLMDQMModule.h:139
void event() override
This method is called for each event.
const KLMElementNumbers * m_ElementNumbers
KLM element numbers.
Definition: KLMDQMModule.h:237
KLMDQMModule()
Constructor.
Definition: KLMDQMModule.cc:23
double m_Revo9DCArrivalTimeMin
Min time for revo9DCArrivalTime for RPC.
Definition: KLMDQMModule.h:136
double m_EKLMScintTimeMax
Max time for EKLM Scint.
Definition: KLMDQMModule.h:133
const int m_ChannelHitHistogramsBKLM
Number of channel hit histograms per sector for BKLM.
Definition: KLMDQMModule.h:112
TH1F * m_TriggerBitsEKLM
Trigger bits: EKLM scintillators.
Definition: KLMDQMModule.h:198
void endRun() override
This method is called if the current run ends.
const EKLMElementNumbers * m_eklmElementNumbers
Element numbers.
Definition: KLMDQMModule.h:240
double m_BKLMScintTimeMax
Max time for BKLM Scint.
Definition: KLMDQMModule.h:127
StoreObjPtr< TRGSummary > m_trgSummary
Trigger summary (event-level L1 bits).
Definition: KLMDQMModule.h:252
void terminate() override
This method is called at the end of the event processing.
TH1F * m_TimeRevo9DCArrivalTime
Time: revo9DCArrivalTime for RPC.
Definition: KLMDQMModule.h:151
TH1F * m_PlaneBKLMZ
Plane occupancy: BKLM, z readout.
Definition: KLMDQMModule.h:163
TH1F * m_TimeScintillatorBKLM
Time: BKLM scintillators.
Definition: KLMDQMModule.h:154
TH1F * m_FE_BKLM_Layer_1
feature extraction status for BKLM Scintillator
Definition: KLMDQMModule.h:222
TH1F * m_DigitsRPC
Number of digits: BKLM RPCs.
Definition: KLMDQMModule.h:180
TH1F * m_TriggersHERInj
Histogram to be used for normalization of occupancy after HER injection.
Definition: KLMDQMModule.h:213
void beginRun() override
Called when entering a new run.
TH1F * m_DigitsAfterHERInj
Number of KLM Digits after LER injection.
Definition: KLMDQMModule.h:210
TH1F * m_DigitsScintillatorBKLM
Number of digits: BKLM scintillators.
Definition: KLMDQMModule.h:183
TH1F * m_TriggersLERInj
Histogram to be used for normalization of occupancy after LER injection.
Definition: KLMDQMModule.h:207
TH1F * m_DigitsMultiStripBKLM
Number of multi-strip digits: BKLM scintillators.
Definition: KLMDQMModule.h:189
const KLMSectorArrayIndex * m_SectorArrayIndex
KLM sector array index.
Definition: KLMDQMModule.h:234
StoreArray< BKLMHit1d > m_BklmHit1ds
BKLM 1d hits.
Definition: KLMDQMModule.h:258
TH1F * m_DigitsAfterLERInj
Number of KLM Digits after LER injection.
Definition: KLMDQMModule.h:204
TH1F * m_DAQInclusion
KLM DAQ inclusion.
Definition: KLMDQMModule.h:145
TH1F * m_FE_EKLM_Plane_0
Standard Readout status for EKLM.
Definition: KLMDQMModule.h:225
StoreArray< RawKLM > m_RawKlms
Raw KLM.
Definition: KLMDQMModule.h:249
KLMTime * m_klmTime
KLM time conversion (for revo9DCArrivalTime).
Definition: KLMDQMModule.h:243
TH1F * m_PlaneEKLM
Plane occupancy: EKLM.
Definition: KLMDQMModule.h:166
TH2F ** m_Spatial2DHitsEKLM[EKLMElementNumbers::getMaximalSectionNumber()]
Spatial distribution of EKLM 2d hits per layer.
Definition: KLMDQMModule.h:216
TH1F * m_EventBackgroundTriggerSummary
Event-level background trigger summary (from TRGSummary).
Definition: KLMDQMModule.h:201
TH1F * m_MaskedChannelsPerSector
Masked channels per sector.
Definition: KLMDQMModule.h:174
double m_EKLMScintTimeMin
Min time for EKLM Scint.
Definition: KLMDQMModule.h:130
const KLMChannelArrayIndex * m_ChannelArrayIndex
KLM channel array index.
Definition: KLMDQMModule.h:231
static constexpr std::array< TRGSummary::ETimingType, 3 > c_KlmL1Triggers
L1 timing trigger bits of interest for KLM DQM (event-level).
Definition: KLMDQMModule.h:264
std::string m_HistogramDirectoryName
Directory for KLM DQM histograms in ROOT file.
Definition: KLMDQMModule.h:142
TH1F * m_FE_BKLM_Layer_0
Standard Readout status for BKLM Scintillator.
Definition: KLMDQMModule.h:219
void defineHisto() override
Definition of the histograms.
Definition: KLMDQMModule.cc:99
StoreArray< KLMHit2d > m_Hit2ds
KLM 2d hits.
Definition: KLMDQMModule.h:261
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.
KLM time conversion.
Definition: KLMTime.h:27
double getCTimePeriod() const
Get CTIME period.
Definition: KLMTime.h:53
void updateConstants()
Update constants from database objects.
Definition: KLMTime.cc:20
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.