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