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