Belle II Software  release-05-01-25
DQMHistAnalysisKLM.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, Leo Piilonen, Vipin Gaur, *
7  * Giacomo De Pietro *
8  * *
9  * This software is provided "as is" without any warranty. *
10  **************************************************************************/
11 
12 /* Own header. */
13 #include <dqm/analysis/modules/DQMHistAnalysisKLM.h>
14 
15 /* Belle 2 headers. */
16 #include <klm/dataobjects/KLMChannelIndex.h>
17 
18 /* ROOT headers. */
19 #include <TClass.h>
20 #include <TROOT.h>
21 
22 /* C++ headers. */
23 #include <algorithm>
24 
25 using namespace Belle2;
26 
27 REG_MODULE(DQMHistAnalysisKLM)
28 
31  m_ProcessedEvents{0},
32  m_ChannelArrayIndex{&(KLMChannelArrayIndex::Instance())},
33  m_SectorArrayIndex{&(KLMSectorArrayIndex::Instance())},
34  m_ElementNumbers{&(KLMElementNumbers::Instance())},
35  m_EklmElementNumbers{&(EKLMElementNumbers::Instance())}
36 {
37  setDescription("Module used to analyze KLM DQM histograms.");
38  addParam("ThresholdForMasked", m_ThresholdForMasked,
39  "Threshold X for masked channels: if a channel has an occupancy X times larger than the average, it will be masked.", 100);
40  addParam("ThresholdForHot", m_ThresholdForHot,
41  "Threshold Y for hot channels: if a channel has an occupancy Y times larger than the average, it will be marked as hot (but not masked).",
42  10);
43  addParam("MinHitsForFlagging", m_MinHitsForFlagging, "Minimal number of hits in a channel required to flag it as 'Masked' or 'Hot'",
44  50);
45  addParam("MinProcessedEventsForMessages", m_MinProcessedEventsForMessagesInput,
46  "Minimal number of processed events required to print error messages", 10000.);
47 
48  m_MinProcessedEventsForMessages = m_MinProcessedEventsForMessagesInput;
49  m_PlaneLine.SetLineColor(kMagenta);
50  m_PlaneLine.SetLineWidth(1);
51  m_PlaneLine.SetLineStyle(2); // dashed
52  m_PlaneText.SetTextAlign(22); // centered, middle
53  m_PlaneText.SetTextColor(kMagenta);
54  m_PlaneText.SetTextFont(42); // Helvetica regular
55  m_PlaneText.SetTextSize(0.02); // 2% of TPad's full height
56 }
57 
59 {
60 }
61 
63 {
65  B2FATAL("The threshold used for hot channels is larger than the one for masked channels."
66  << LogVar("Threshold for hot channels", m_ThresholdForHot)
67  << LogVar("Threshold for masked channels", m_ThresholdForMasked));
68 }
69 
71 {
72 }
73 
75 {
76  if (!m_ElectronicsMap.isValid())
77  B2FATAL("No KLM electronics map.");
79  m_ProcessedEvents = 0.;
80  m_DeadBarrelModules.clear();
81  m_DeadEndcapModules.clear();
82  m_MaskedChannels.clear();
83 }
84 
86 {
87 }
88 
90 {
91  TH1* histogram = findHist("DAQ/Nevent");
92  if (histogram == nullptr) {
93  B2WARNING("DAQ DQM histogram DAQ/Nevent is not found.");
94  /* Set the minimal number of processed events to 0 if we can't determine the processed events. */
96  return 0.;
97  }
98  return histogram->GetEntries();
99 }
100 
102  int subdetector, int section, int sector,
103  TH1* histogram, TCanvas* canvas, TLatex& latex)
104 {
105  double x = 0.15;
106  double y = 0.85;
107  int i, n;
108  std::map<uint16_t, double> moduleHitMap;
109  std::map<uint16_t, double>::iterator it;
110  double average = 0;
111  int channelSubdetector, channelSection, channelSector;
112  int layer, plane, strip;
113  std::string str;
114  canvas->Clear();
115  canvas->cd();
116  histogram->SetStats(false);
117  histogram->Draw();
118  n = histogram->GetXaxis()->GetNbins();
119  for (i = 1; i <= n; i++) {
120  uint16_t channelIndex = std::round(histogram->GetBinCenter(i));
121  uint16_t channelNumber = m_ChannelArrayIndex->getNumber(channelIndex);
122  double nHitsPerModule = histogram->GetBinContent(i);
123  average = average + nHitsPerModule;
125  channelNumber, &channelSubdetector, &channelSection, &channelSector,
126  &layer, &plane, &strip);
127  if ((channelSubdetector != subdetector) ||
128  (channelSection != section) ||
129  (channelSector != sector))
130  B2FATAL("Inconsistent element numbers.");
131  uint16_t module = m_ElementNumbers->moduleNumber(
132  subdetector, section, sector, layer);
133  it = moduleHitMap.find(module);
134  if (it == moduleHitMap.end())
135  moduleHitMap.insert(std::pair<uint16_t, double>(module, nHitsPerModule));
136  else
137  it->second += nHitsPerModule;
138  }
139  unsigned int activeModuleChannels = 0;
140  for (it = moduleHitMap.begin(); it != moduleHitMap.end(); ++it) {
141  uint16_t moduleNumber = it->first;
142  if (it->second != 0) {
143  activeModuleChannels += m_ElementNumbers->getNChannelsModule(moduleNumber);
144  continue;
145  }
147  moduleNumber, &channelSubdetector, &channelSection, &channelSector, &layer);
148  /* Channel with plane = 1, strip = 1 exists for any BKLM or EKLM module. */
149  uint16_t channel = m_ElementNumbers->channelNumber(
150  channelSubdetector, channelSection, channelSector,
151  layer, 1, 1);
152  const KLMElectronicsChannel* electronicsChannel =
153  m_ElectronicsMap->getElectronicsChannel(channel);
154  if (electronicsChannel == nullptr)
155  B2FATAL("Incomplete KLM electronics map.");
156  str = "No data from HSLB ";
157  if (channelSubdetector == KLMElementNumbers::c_BKLM) {
158  str += BKLMElementNumbers::getHSLBName(electronicsChannel->getCopper(),
159  electronicsChannel->getSlot());
160  } else {
161  str += EKLMElementNumbers::getHSLBName(electronicsChannel->getCopper(),
162  electronicsChannel->getSlot());
163  }
164  str += ", lane " + std::to_string(electronicsChannel->getLane());
165  latex.DrawLatexNDC(x, y, str.c_str());
166  y -= 0.05;
167  /* Store the module number, used later in processPlaneHistogram
168  * to color the canvas with red and to raise up an alarm. */
169  if (channelSubdetector == KLMElementNumbers::c_BKLM) {
170  std::vector<uint16_t>::iterator ite = std::find(m_DeadBarrelModules.begin(),
171  m_DeadBarrelModules.end(),
172  moduleNumber);
173  if (ite == m_DeadBarrelModules.end())
174  m_DeadBarrelModules.push_back(moduleNumber);
175  } else {
176  std::vector<uint16_t>::iterator ite = std::find(m_DeadEndcapModules.begin(),
177  m_DeadEndcapModules.end(),
178  moduleNumber);
179  if (ite == m_DeadEndcapModules.end())
180  m_DeadEndcapModules.push_back(moduleNumber);
181  }
182  }
183  if (activeModuleChannels == 0)
184  return;
185  average /= activeModuleChannels;
186  for (i = 1; i <= n; ++i) {
187  uint16_t channelIndex = std::round(histogram->GetBinCenter(i));
188  uint16_t channelNumber = m_ChannelArrayIndex->getNumber(channelIndex);
189  double nHits = histogram->GetBinContent(i);
191  channelNumber, &channelSubdetector, &channelSection, &channelSector,
192  &layer, &plane, &strip);
193  std::string channelStatus = "Normal";
194  if ((nHits > average * m_ThresholdForMasked) && (nHits > m_MinHitsForFlagging)) {
195  channelStatus = "Masked";
196  std::vector<uint16_t>::iterator ite = std::find(m_MaskedChannels.begin(),
197  m_MaskedChannels.end(),
198  channelNumber);
199  if (ite == m_MaskedChannels.end())
200  m_MaskedChannels.push_back(channelNumber);
201  B2DEBUG(20, "KLM@MaskMe " << channelNumber);
202  } else if ((nHits > average * m_ThresholdForHot) && (nHits > m_MinHitsForFlagging)) {
203  channelStatus = "Hot";
204  }
205  if (channelStatus != "Normal") {
206  const KLMElectronicsChannel* electronicsChannel =
207  m_ElectronicsMap->getElectronicsChannel(channelNumber);
208  if (electronicsChannel == nullptr)
209  B2FATAL("Incomplete BKLM electronics map.");
210  if (channelStatus == "Masked")
211  histogram->SetBinContent(i, 0);
212  str = channelStatus + " channel: HSLB ";
213  if (channelSubdetector == KLMElementNumbers::c_BKLM) {
214  str += BKLMElementNumbers::getHSLBName(electronicsChannel->getCopper(),
215  electronicsChannel->getSlot());
216  } else {
217  str += EKLMElementNumbers::getHSLBName(electronicsChannel->getCopper(),
218  electronicsChannel->getSlot());
219  }
220  str += (", lane " + std::to_string(electronicsChannel->getLane()) +
221  ", axis " + std::to_string(electronicsChannel->getAxis()) +
222  ", channel " + std::to_string(electronicsChannel->getChannel()));
223  latex.DrawLatexNDC(x, y, str.c_str());
224  y -= 0.05;
225  }
226  }
227  canvas->Modified();
228 }
229 
231  const std::string& histName)
232 {
233  TH1* histogram = findHist("KLM/" + histName);
234  if (histogram == nullptr) {
235  B2ERROR("KLM DQM histogram KLM/" << histName << " is not found.");
236  return;
237  }
238  TCanvas* canvas = findCanvas("KLM/c_" + histName);
239  if (canvas == nullptr) {
240  B2ERROR("KLM DQM histogram canvas KLM/c_" << histName << " is not found.");
241  return;
242  }
243  histogram->Clear();
244  canvas->Clear();
245  canvas->cd();
246  if (m_MaskedChannels.size() > 0) {
247  int channelSubdetector, channelSection, channelSector;
248  int layer, plane, strip;
249  for (uint16_t channel : m_MaskedChannels) {
251  channel, &channelSubdetector, &channelSection, &channelSector,
252  &layer, &plane, &strip);
253  uint16_t sectorNumber;
254  if (channelSubdetector == KLMElementNumbers::c_BKLM)
255  sectorNumber = m_ElementNumbers->sectorNumberBKLM(channelSection, channelSector);
256  else
257  sectorNumber = m_ElementNumbers->sectorNumberEKLM(channelSection, channelSector);
258  uint16_t sectorIndex = m_SectorArrayIndex->getIndex(sectorNumber);
259  histogram->Fill(sectorIndex);
260  }
261  }
262  histogram->SetStats(false);
263  histogram->Draw();
264  canvas->Modified();
265 }
266 
268  const std::string& histName, TLatex& latex)
269 {
270  std::string name, alarm;
271  const double histMinNDC = 0.1;
272  const double histMaxNDC = 0.9;
273  const double histRangeNDC = histMaxNDC - histMinNDC;
274  int moduleSubdetector, moduleSection, moduleSector, moduleLayer;
275  double xAlarm = 0.15;
276  double yAlarm = 0.8;
277  TH1* histogram = findHist("KLM/" + histName);
278  if (histogram == nullptr) {
279  B2ERROR("KLM DQM histogram KLM/" << histName << " is not found.");
280  return;
281  }
282  TCanvas* canvas = findCanvas("KLM/c_" + histName);
283  if (canvas == nullptr) {
284  B2ERROR("KLM DQM histogram canvas KLM/c_" << histName << " is not found.");
285  return;
286  }
287  canvas->Clear();
288  canvas->cd();
289  histogram->SetStats(false);
290  histogram->Draw();
291  if (histName.find("bklm") != std::string::npos) {
292  /* First draw the vertical lines and the sector names. */
293  const double maximalSector = BKLMElementNumbers::getMaximalSectorGlobalNumber();
294  for (int sector = 0; sector < BKLMElementNumbers::getMaximalSectorGlobalNumber(); ++sector) {
295  double xLineNDC = histMinNDC + (histRangeNDC * sector) / maximalSector;
296  double xTextNDC = histMinNDC + (histRangeNDC * (sector + 0.5)) / maximalSector;
297  double yTextNDC = histMinNDC + 0.98 * histRangeNDC;
298  if (sector > 0)
299  m_PlaneLine.DrawLineNDC(xLineNDC, histMinNDC, xLineNDC, histMaxNDC);
300  name = "B";
301  if (sector < 8)
302  name += "B";
303  else
304  name += "F";
305  name += std::to_string(sector % 8);
306  m_PlaneText.DrawTextNDC(xTextNDC, yTextNDC, name.c_str());
307  }
308  /* Then, color the canvas with red if there is a dead module
309  * and write an error message. */
310  if (m_DeadBarrelModules.size() == 0) {
311  canvas->Pad()->SetFillColor(kWhite);
313  for (uint16_t module : m_DeadBarrelModules) {
315  module, &moduleSubdetector, &moduleSection, &moduleSector, &moduleLayer);
316  alarm = "No data from " + m_ElementNumbers->getSectorDAQName(moduleSubdetector, moduleSection, moduleSector);
317  alarm += ", layer " + std::to_string(moduleLayer);
318  latex.DrawLatexNDC(xAlarm, yAlarm, alarm.c_str());
319  yAlarm -= 0.05;
320  }
321  alarm = "Call the KLM experts immediately!";
322  latex.DrawLatexNDC(xAlarm, yAlarm, alarm.c_str());
323  canvas->Pad()->SetFillColor(kRed);
324  }
325  } else {
326  /* First draw the vertical lines and the sector names. */
327  const double maximalLayer = EKLMElementNumbers::getMaximalLayerGlobalNumber();
328  for (int layerGlobal = 1; layerGlobal <= maximalLayer; ++layerGlobal) {
329  double xLineNDC = histMinNDC + (histRangeNDC * layerGlobal) / maximalLayer;
330  double xTextNDC = histMinNDC + (histRangeNDC * (layerGlobal - 0.5)) / maximalLayer;
331  double yTextNDC = histMinNDC + 0.98 * histRangeNDC;
332  if (layerGlobal < maximalLayer)
333  m_PlaneLine.DrawLineNDC(xLineNDC, histMinNDC, xLineNDC, histMaxNDC);
334  int section, layer;
336  layerGlobal, &section, &layer);
338  name = "B";
339  else
340  name = "F";
341  name += std::to_string(layer);
342  m_PlaneText.DrawTextNDC(xTextNDC, yTextNDC, name.c_str());
343  }
344  /* Then, color the canvas with red if there is a dead module
345  * and write an error message. */
346  if (m_DeadEndcapModules.size() == 0) {
347  canvas->Pad()->SetFillColor(kWhite);
349  for (uint16_t module : m_DeadEndcapModules) {
351  module, &moduleSubdetector, &moduleSection, &moduleSector, &moduleLayer);
352  alarm = "No data from " + m_ElementNumbers->getSectorDAQName(moduleSubdetector, moduleSection, moduleSector);
353  alarm += ", layer " + std::to_string(moduleLayer);
354  latex.DrawLatexNDC(xAlarm, yAlarm, alarm.c_str());
355  yAlarm -= 0.05;
356  }
357  alarm = "Call the KLM experts immediately!";
358  latex.DrawLatexNDC(xAlarm, yAlarm, alarm.c_str());
359  canvas->Pad()->SetFillColor(kRed);
360  }
361  }
362  canvas->Modified();
363  canvas->Update();
364 }
365 
366 TCanvas* DQMHistAnalysisKLMModule::findCanvas(const std::string& canvasName)
367 {
368  TIter nextkey(gROOT->GetListOfCanvases());
369  TObject* obj = nullptr;
370  while ((obj = dynamic_cast<TObject*>(nextkey()))) {
371  if (obj->IsA()->InheritsFrom("TCanvas")) {
372  if (obj->GetName() == canvasName)
373  return dynamic_cast<TCanvas*>(obj);
374  }
375  }
376  return nullptr;
377 }
378 
379 
381 {
382  int isKlmIncluded = 1;
383  /* If KLM is not included, stop here and return. */
384  TH1* daqInclusion = findHist("KLM/daq_inclusion");
385  if (not(daqInclusion == nullptr)) {
386  isKlmIncluded = daqInclusion->GetBinContent(daqInclusion->GetXaxis()->FindBin("Yes"));
387  if (isKlmIncluded == 0)
388  return;
389  }
390  /* Make sure that the vectors are cleared at each DQM refresh. */
391  m_DeadBarrelModules.clear();
392  m_DeadEndcapModules.clear();
393  m_MaskedChannels.clear();
395  std::string str, histogramName, canvasName;
396  TLatex latex;
397  latex.SetTextColor(kRed);
398  latex.SetTextAlign(11);
400  for (KLMChannelIndex& klmSector : klmSectors) {
401  int nHistograms;
402  if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM)
403  nHistograms = 2;
404  else
405  nHistograms = 3;
406  for (int j = 0; j < nHistograms; j++) {
407  str = "strip_hits_subdetector_" +
408  std::to_string(klmSector.getSubdetector()) +
409  "_section_" + std::to_string(klmSector.getSection()) +
410  "_sector_" + std::to_string(klmSector.getSector()) +
411  "_" + std::to_string(j);
412  histogramName = "KLM/" + str;
413  canvasName = "KLM/c_" + str;
414  TH1* histogram = findHist(histogramName);
415  if (histogram == nullptr) {
416  B2ERROR("KLM DQM histogram " << histogramName << " is not found.");
417  continue;
418  }
419  TCanvas* canvas = findCanvas(canvasName);
420  if (canvas == nullptr) {
421  B2ERROR("KLM DQM histogram canvas " << canvasName << " is not found.");
422  continue;
423  }
425  klmSector.getSubdetector(), klmSector.getSection(),
426  klmSector.getSector(), histogram, canvas, latex);
427  }
428  }
429  fillMaskedChannelsHistogram("masked_channels");
430  latex.SetTextColor(kBlue);
431  processPlaneHistogram("plane_bklm_phi", latex);
432  processPlaneHistogram("plane_bklm_z", latex);
433  processPlaneHistogram("plane_eklm", latex);
434 }
Belle2::KLMElementNumbers::sectorNumberBKLM
uint16_t sectorNumberBKLM(int section, int sector) const
Get sector number for BKLM.
Definition: KLMElementNumbers.cc:222
Belle2::KLMElectronicsChannel::getChannel
int getChannel() const
Get channel.
Definition: KLMElectronicsChannel.h:145
Belle2::DQMHistAnalysisKLMModule::m_ThresholdForMasked
int m_ThresholdForMasked
Threshold for masked channels.
Definition: DQMHistAnalysisKLM.h:135
Belle2::DQMHistAnalysisKLMModule::m_ChannelArrayIndex
const KLMChannelArrayIndex * m_ChannelArrayIndex
KLM channel array index.
Definition: DQMHistAnalysisKLM.h:165
Belle2::DQMHistAnalysisKLMModule::m_MaskedChannels
std::vector< uint16_t > m_MaskedChannels
Vector of masked channels.
Definition: DQMHistAnalysisKLM.h:156
Belle2::DQMHistAnalysisKLMModule::processPlaneHistogram
void processPlaneHistogram(const std::string &histName, TLatex &latex)
Process histogram containing the number of hits in plane.
Definition: DQMHistAnalysisKLM.cc:267
Belle2::DQMHistAnalysisKLMModule::~DQMHistAnalysisKLMModule
~DQMHistAnalysisKLMModule()
Destructor.
Definition: DQMHistAnalysisKLM.cc:58
Belle2::EKLMElementNumbers::getHSLBName
static std::string getHSLBName(int copper, int slot)
Get HSLB name.
Definition: EKLMElementNumbers.cc:273
Belle2::DQMHistAnalysisKLMModule::m_PlaneLine
TLine m_PlaneLine
TLine for boundary in plane histograms.
Definition: DQMHistAnalysisKLM.h:159
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::BKLMElementNumbers::getHSLBName
static std::string getHSLBName(int copper, int slot)
Get HSLB name.
Definition: BKLMElementNumbers.cc:203
Belle2::DQMHistAnalysisKLMModule::m_ElectronicsMap
DBObjPtr< KLMElectronicsMap > m_ElectronicsMap
Electronics map.
Definition: DQMHistAnalysisKLM.h:177
Belle2::KLMElectronicsChannel::getSlot
int getSlot() const
Get slot.
Definition: KLMElectronicsChannel.h:94
Belle2::DQMHistAnalysisKLMModule::m_MinHitsForFlagging
int m_MinHitsForFlagging
Minimal number of hits for flagging.
Definition: DQMHistAnalysisKLM.h:141
Belle2::KLMSectorArrayIndex::Instance
static const KLMSectorArrayIndex & Instance()
Instantiation.
Definition: KLMSectorArrayIndex.cc:28
Belle2::DQMHistAnalysisKLMModule::event
void event() override
This method is called for each event.
Definition: DQMHistAnalysisKLM.cc:380
Belle2::KLMElementNumbers::channelNumber
uint16_t channelNumber(int subdetector, int section, int sector, int layer, int plane, int strip) const
Get channel number.
Definition: KLMElementNumbers.cc:37
Belle2::KLMElementNumbers::Instance
static const KLMElementNumbers & Instance()
Instantiation.
Definition: KLMElementNumbers.cc:31
Belle2::DQMHistAnalysisKLMModule::m_MinProcessedEventsForMessages
double m_MinProcessedEventsForMessages
Minimal number of processed events for error messages.
Definition: DQMHistAnalysisKLM.h:147
Belle2::KLMElementNumbers::moduleNumberToElementNumbers
void moduleNumberToElementNumbers(uint16_t module, int *subdetector, int *section, int *sector, int *layer) const
Get element numbers by module number.
Definition: KLMElementNumbers.cc:186
Belle2::DQMHistAnalysisModule::findHist
static TH1 * findHist(const std::string &histname)
Find histogram.
Definition: DQMHistAnalysis.cc:83
Belle2::EKLMElementNumbers::c_BackwardSection
@ c_BackwardSection
Backward.
Definition: EKLMElementNumbers.h:44
Belle2::EKLMElementNumbers::getMaximalLayerGlobalNumber
static constexpr int getMaximalLayerGlobalNumber()
Get maximal detector layer global number.
Definition: EKLMElementNumbers.h:337
Belle2::KLMElectronicsChannel::getAxis
int getAxis() const
Get axis.
Definition: KLMElectronicsChannel.h:128
Belle2::KLMChannelIndex::c_IndexLevelSector
@ c_IndexLevelSector
Sector.
Definition: KLMChannelIndex.h:49
Belle2::KLMElementArrayIndex::getNumber
uint16_t getNumber(uint16_t index) const
Get element number.
Definition: KLMElementArrayIndex.cc:63
Belle2::DQMHistAnalysisKLMModule::getProcessedEvents
double getProcessedEvents()
Get number of processed events.
Definition: DQMHistAnalysisKLM.cc:89
Belle2::moduleHitMap
TH2 * moduleHitMap(TH1 *hitMap, int moduleID)
Make hit map in HAPD view (12*12 channels)
Definition: hitMapMaker.cc:42
Belle2::DQMHistAnalysisKLMModule::m_MinProcessedEventsForMessagesInput
double m_MinProcessedEventsForMessagesInput
Input parameter for minimal number of processed events for error messages.
Definition: DQMHistAnalysisKLM.h:144
Belle2::DQMHistAnalysisKLMModule::m_EklmElementNumbers
const EKLMElementNumbers * m_EklmElementNumbers
EKLM element numbers.
Definition: DQMHistAnalysisKLM.h:174
Belle2::KLMElementNumbers::getNChannelsModule
unsigned int getNChannelsModule(uint16_t module) const
Get number of channels in module.
Definition: KLMElementNumbers.cc:208
Belle2::EKLMElementNumbers::layerNumberToElementNumbers
void layerNumberToElementNumbers(int layerGlobal, int *section, int *layer) const
Get element numbers by detector layer global number.
Definition: EKLMElementNumbers.cc:130
Belle2::KLMElementNumbers::sectorNumberEKLM
uint16_t sectorNumberEKLM(int section, int sector) const
Get sector number for EKLM.
Definition: KLMElementNumbers.cc:229
Belle2::DQMHistAnalysisKLMModule::m_ThresholdForHot
int m_ThresholdForHot
Threshold for hot channels.
Definition: DQMHistAnalysisKLM.h:138
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisKLMModule::endRun
void endRun() override
This method is called if the current run ends.
Definition: DQMHistAnalysisKLM.cc:85
Belle2::BKLMElementNumbers::getMaximalSectorGlobalNumber
static constexpr int getMaximalSectorGlobalNumber()
Get maximal sector global number.
Definition: BKLMElementNumbers.h:267
Belle2::DQMHistAnalysisKLMModule::initialize
void initialize() override
Initializer.
Definition: DQMHistAnalysisKLM.cc:62
Belle2::DQMHistAnalysisKLMModule::m_DeadEndcapModules
std::vector< uint16_t > m_DeadEndcapModules
Vector of dead endcap modules.
Definition: DQMHistAnalysisKLM.h:153
LogVar
Class to store variables with their name which were sent to the logging service.
Definition: LogVariableStream.h:24
Belle2::EKLMElementNumbers::Instance
static const EKLMElementNumbers & Instance()
Instantiation.
Definition: EKLMElementNumbers.cc:20
Belle2::DQMHistAnalysisKLMModule::findCanvas
TCanvas * findCanvas(const std::string &canvasName)
Find TCanvas that matches a given name.
Definition: DQMHistAnalysisKLM.cc:366
Belle2::KLMElementNumbers::c_BKLM
@ c_BKLM
BKLM.
Definition: KLMElementNumbers.h:47
Belle2::DQMHistAnalysisKLMModule
Analysis of KLM DQM histograms.
Definition: DQMHistAnalysisKLM.h:52
Belle2::DQMHistAnalysisKLMModule::analyseChannelHitHistogram
void analyseChannelHitHistogram(int subdetector, int section, int sector, TH1 *histogram, TCanvas *canvas, TLatex &latex)
Analyse channel hit histogram.
Definition: DQMHistAnalysisKLM.cc:101
Belle2::DQMHistAnalysisKLMModule::terminate
void terminate() override
This method is called at the end of the event processing.
Definition: DQMHistAnalysisKLM.cc:70
Belle2::DQMHistAnalysisKLMModule::m_ElementNumbers
const KLMElementNumbers * m_ElementNumbers
KLM element numbers.
Definition: DQMHistAnalysisKLM.h:171
Belle2::KLMChannelArrayIndex::Instance
static const KLMChannelArrayIndex & Instance()
Instantiation.
Definition: KLMChannelArrayIndex.cc:28
Belle2::KLMElectronicsChannel::getCopper
int getCopper() const
Get copper.
Definition: KLMElectronicsChannel.h:77
Belle2::KLMChannelIndex
KLM channel index.
Definition: KLMChannelIndex.h:33
Belle2::KLMElectronicsChannel
BKLM electronics channel.
Definition: KLMElectronicsChannel.h:33
Belle2::DQMHistAnalysisKLMModule::m_SectorArrayIndex
const KLMSectorArrayIndex * m_SectorArrayIndex
KLM sector array index.
Definition: DQMHistAnalysisKLM.h:168
Belle2::DQMHistAnalysisKLMModule::m_DeadBarrelModules
std::vector< uint16_t > m_DeadBarrelModules
Vector of dead barrel modules.
Definition: DQMHistAnalysisKLM.h:150
Belle2::DQMHistAnalysisKLMModule::m_PlaneText
TText m_PlaneText
TText for names in plane histograms.
Definition: DQMHistAnalysisKLM.h:162
Belle2::KLMElectronicsChannel::getLane
int getLane() const
Get lane.
Definition: KLMElectronicsChannel.h:111
Belle2::KLMElementNumbers::moduleNumber
uint16_t moduleNumber(int subdetector, int section, int sector, int layer) const
Get module number.
Definition: KLMElementNumbers.cc:149
Belle2::DQMHistAnalysisKLMModule::beginRun
void beginRun() override
Called when entering a new run.
Definition: DQMHistAnalysisKLM.cc:74
Belle2::DQMHistAnalysisKLMModule::fillMaskedChannelsHistogram
void fillMaskedChannelsHistogram(const std::string &histName)
Fill histogram containing masked channels per sector.
Definition: DQMHistAnalysisKLM.cc:230
Belle2::DQMHistAnalysisKLMModule::m_ProcessedEvents
double m_ProcessedEvents
Number of processed events.
Definition: DQMHistAnalysisKLM.h:132
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::KLMElementArrayIndex::getIndex
uint16_t getIndex(uint16_t number) const
Get element index.
Definition: KLMElementArrayIndex.cc:54
Belle2::KLMElementNumbers::channelNumberToElementNumbers
void channelNumberToElementNumbers(uint16_t channel, int *subdetector, int *section, int *sector, int *layer, int *plane, int *strip) const
Get element numbers by channel number.
Definition: KLMElementNumbers.cc:106
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27