Belle II Software development
DQMHistAnalysisKLM.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 <dqm/analysis/modules/DQMHistAnalysisKLM.h>
11
12/* Basf2 headers. */
13#include <klm/dataobjects/KLMChannelIndex.h>
14
15/* ROOT headers. */
16#include <TClass.h>
17#include <TROOT.h>
18#include <TStyle.h>
19
20/* C++ headers. */
21#include <algorithm>
22
23using namespace Belle2;
24
25REG_MODULE(DQMHistAnalysisKLM);
26
29 m_ProcessedEvents{0},
30 m_IsNullRun{false},
31 m_ChannelArrayIndex{&(KLMChannelArrayIndex::Instance())},
32 m_SectorArrayIndex{&(KLMSectorArrayIndex::Instance())},
33 m_ElementNumbers{&(KLMElementNumbers::Instance())},
34 m_EklmElementNumbers{&(EKLMElementNumbers::Instance())}
35{
36 setDescription("Module used to analyze KLM DQM histograms.");
37 addParam("ThresholdForMasked", m_ThresholdForMasked,
38 "Threshold X for masked channels: if a channel has an occupancy X times larger than the average, it will be masked.", 100);
39 addParam("ThresholdForHot", m_ThresholdForHot,
40 "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).",
41 10);
42 addParam("ThresholdForLog", m_ThresholdForLog,
43 "Threshold Z for log scale view: if a channel has an occupancy Z times larger than the average, canvas shifts to log scale.",
44 20);
45 addParam("MinHitsForFlagging", m_MinHitsForFlagging, "Minimal number of hits in a channel required to flag it as 'Masked' or 'Hot'",
46 50);
47 addParam("MinProcessedEventsForMessages", m_MinProcessedEventsForMessagesInput,
48 "Minimal number of processed events required to print error messages", 10000.);
49 addParam("MinEntries", m_minEntries,
50 "Minimal number for delta histogram updates", 50000.);
51 addParam("MessageThreshold", m_MessageThreshold,
52 "Max number of messages to show up in channel occupancy plots", 12);
53 addParam("HistogramDirectoryName", m_histogramDirectoryName, "Name of histogram directory", std::string("KLM"));
54
56 m_2DHitsLine.SetLineColor(kRed);
57 m_2DHitsLine.SetLineWidth(3);
58 m_2DHitsLine.SetLineStyle(2); // dashed
59 m_PlaneLine.SetLineColor(kMagenta);
60 m_PlaneLine.SetLineWidth(1);
61 m_PlaneLine.SetLineStyle(2); // dashed
62 m_PlaneText.SetTextAlign(22); // centered, middle
63 m_PlaneText.SetTextColor(kMagenta);
64 m_PlaneText.SetTextFont(42); // Helvetica regular
65 m_PlaneText.SetTextSize(0.02); // 2% of TPad's full height
66}
67
68
70{
72
74 B2FATAL("The threshold used for hot channels is larger than the one for masked channels."
75 << LogVar("Threshold for hot channels", m_ThresholdForHot)
76 << LogVar("Threshold for masked channels", m_ThresholdForMasked));
77
78 // register plots for delta histogramming
79 addDeltaPar(m_histogramDirectoryName, "time_rpc", HistDelta::c_Entries, m_minEntries, 1);
80 addDeltaPar(m_histogramDirectoryName, "time_scintillator_bklm", HistDelta::c_Entries, m_minEntries, 1);
81 addDeltaPar(m_histogramDirectoryName, "time_scintillator_eklm", HistDelta::c_Entries, m_minEntries, 1);
82
83 //register EPICS PVs
84 registerEpicsPV("KLM:MaskedChannels", "MaskedChannels");
85 registerEpicsPV("KLM:DeadBarrelModules", "DeadBarrelModules");
86 registerEpicsPV("KLM:DeadEndcapModules", "DeadEndcapModules");
87
88 gROOT->cd();
89
90 std::string str;
92 for (KLMChannelIndex& klmSector : klmIndex) {
93 int nHistograms;
94 if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM)
95 nHistograms = 2;
96 else
97 nHistograms = 3;
98 for (int j = 0; j < nHistograms; j++) {
99 str = "strip_hits_subdetector_" +
100 std::to_string(klmSector.getSubdetector()) +
101 "_section_" + std::to_string(klmSector.getSection()) +
102 "_sector_" + std::to_string(klmSector.getSector()) +
103 "_" + std::to_string(j);
104 addDeltaPar(m_histogramDirectoryName, str, HistDelta::c_Entries, m_minEntries, 1);
105
106 }
107 }
108
109}
110
112{
113}
114
116{
117 if (!m_ElectronicsMap.isValid())
118 B2FATAL("No KLM electronics map.");
121 m_DeadBarrelModules.clear();
122 m_DeadEndcapModules.clear();
123 m_MaskedChannels.clear();
124
125 m_IsNullRun = (getRunType() == "null");
126}
127
129{
130 int hist_max_bin; double max_position;
131 TH1* time_rpc = findHist(m_histogramDirectoryName + "/time_rpc");
132 if (time_rpc) {
133 hist_max_bin = time_rpc->GetMaximumBin();
134 max_position = time_rpc->GetXaxis()->GetBinCenter(hist_max_bin);
135 m_monObj->setVariable("RPC_Time_Peak", max_position);
136 }
137
138 TH1* time_scint_bklm = findHist(m_histogramDirectoryName + "/time_scintillator_bklm");
139 if (time_scint_bklm) {
140 hist_max_bin = time_scint_bklm->GetMaximumBin();
141 max_position = time_scint_bklm->GetXaxis()->GetBinCenter(hist_max_bin);
142 m_monObj->setVariable("BKLM_Scint_Time_Peak", max_position);
143 }
144
145 TH1* time_scint_eklm = findHist(m_histogramDirectoryName + "/time_scintillator_eklm");
146 if (time_scint_eklm) {
147 hist_max_bin = time_scint_eklm->GetMaximumBin();
148 max_position = time_scint_eklm->GetXaxis()->GetBinCenter(hist_max_bin);
149 m_monObj->setVariable("EKLM_Scint_Time_Peak", max_position);
150 }
151}
152
154{
156 B2WARNING("Either DAQ/Nevent is not found or Nevent = 0.");
157 /* Set the minimal number of processed events to 0 if we can't determine the processed events. */
159 }
161}
162
163void DQMHistAnalysisKLMModule::deltaDrawer(TH1* delta, TH1* histogram, TCanvas* canvas)
164{
165 if (delta != nullptr) {
166 Double_t scale = (Double_t) histogram->Integral(); //want delta and histo to have same norm
167
168 // delta != nullptr should take care of whether update condition is met.
169 delta->SetLineColor(kBlackBody); //random choice of not green or blue
170 delta->SetLineStyle(4);
171 delta->DrawNormalized("SAME", scale); //normalize delta to histo
172 canvas->Modified();
173 canvas->Update();
174 }
175}
176
178 int subdetector, int section, int sector, int index,
179 TH1* histogram, TH1* delta, TCanvas* canvas, TLatex& latex)
180{
181 double x = 0.15;
182 double y = 0.85;
183 int i, n;
184 std::map<KLMModuleNumber, double> moduleHitMap;
185 std::map<KLMModuleNumber, double>::iterator it;
186 double average = 0;
187 int channelSubdetector, channelSection, channelSector;
188 int layer, plane, strip;
189 std::string str;
190 canvas->Clear();
191 canvas->cd();
192 histogram->SetStats(false);
193 canvas->SetLogy(0); //initialize to start without logscale
194 histogram->Draw();
195 deltaDrawer(delta, histogram, canvas); //draw normalized delta on top
196 n = histogram->GetXaxis()->GetNbins();
197
198 /* call reference histograms from base class*/
199 TH1* ref_histogram = findRefHist(histogram->GetName(), false);
200 float ref_average = 0;
201
202 if (ref_histogram != nullptr) {
203 for (i = 1; i <= n; i++) {
204 double nHitsPerModuleRef = ref_histogram->GetBinContent(i);
205 ref_average = ref_average + nHitsPerModuleRef;
206 }
207 }
208
209 for (i = 1; i <= n; i++) {
210 KLMChannelNumber channelIndex = std::round(histogram->GetBinCenter(i));
211 KLMChannelNumber channelNumber =
212 m_ChannelArrayIndex->getNumber(channelIndex);
213 double nHitsPerModule = histogram->GetBinContent(i);
214 average = average + nHitsPerModule;
216 channelNumber, &channelSubdetector, &channelSection, &channelSector,
217 &layer, &plane, &strip);
218 if ((channelSubdetector != subdetector) ||
219 (channelSection != section) ||
220 (channelSector != sector))
221 B2FATAL("Inconsistent element numbers.");
223 subdetector, section, sector, layer);
224 it = moduleHitMap.find(module);
225 if (it == moduleHitMap.end()) {
226 moduleHitMap.insert(std::pair<KLMModuleNumber, double>(
227 module, nHitsPerModule));
228 } else {
229 it->second += nHitsPerModule;
230 }
231 }
232 unsigned int activeModuleChannels = 0;
233 int message_counter = 0;
234 for (it = moduleHitMap.begin(); it != moduleHitMap.end(); ++it) {
235 KLMModuleNumber moduleNumber = it->first;
236 if (it->second != 0) {
237 activeModuleChannels += m_ElementNumbers->getNChannelsModule(moduleNumber);
238 continue;
239 }
241 moduleNumber, &channelSubdetector, &channelSection, &channelSector, &layer);
242 /* Channel with plane = 1, strip = 1 exists for any BKLM or EKLM module. */
243 KLMChannelNumber channel =
245 channelSubdetector, channelSection, channelSector, layer, 1, 1);
246 const KLMElectronicsChannel* electronicsChannel =
247 m_ElectronicsMap->getElectronicsChannel(channel);
248 if (electronicsChannel == nullptr)
249 B2FATAL("Incomplete KLM electronics map.");
250 str = "No data from lane " + std::to_string(electronicsChannel->getLane());
251 latex.DrawLatexNDC(x, y, str.c_str());
252 y -= 0.05;
253 message_counter++;
254 /* Store the module number, used later in processPlaneHistogram
255 * to color the canvas with red and to raise up an alarm. */
256 if (channelSubdetector == KLMElementNumbers::c_BKLM) {
257 std::vector<KLMModuleNumber>::iterator ite =
258 std::find(m_DeadBarrelModules.begin(),
260 moduleNumber);
261 if (ite == m_DeadBarrelModules.end())
262 m_DeadBarrelModules.push_back(moduleNumber);
263 } else {
264 std::vector<KLMModuleNumber>::iterator ite = std::find(m_DeadEndcapModules.begin(),
266 moduleNumber);
267 if (ite == m_DeadEndcapModules.end())
268 m_DeadEndcapModules.push_back(moduleNumber);
269 }
270 }
271 if (activeModuleChannels == 0)
272 return;
273 average /= activeModuleChannels;
274 ref_average /= activeModuleChannels;
275
276 for (i = 1; i <= n; ++i) {
277 KLMChannelNumber channelIndex = std::round(histogram->GetBinCenter(i));
278 KLMChannelNumber channelNumber =
279 m_ChannelArrayIndex->getNumber(channelIndex);
280 double nHits = histogram->GetBinContent(i);
282 channelNumber, &channelSubdetector, &channelSection, &channelSector,
283 &layer, &plane, &strip);
284 std::string channelStatus = "Normal";
285 if ((nHits > average * m_ThresholdForMasked) && (nHits > m_MinHitsForFlagging)) {
286 channelStatus = "Masked";
287 std::vector<KLMModuleNumber>::iterator ite =
288 std::find(m_MaskedChannels.begin(),
289 m_MaskedChannels.end(),
290 channelNumber);
291 if (ite == m_MaskedChannels.end())
292 m_MaskedChannels.push_back(channelNumber);
293 B2DEBUG(20, "KLM@MaskMe " << channelNumber);
294 } else if ((nHits > average * m_ThresholdForHot) && (nHits > m_MinHitsForFlagging)) {
295 channelStatus = "Hot";
296 }
297 if (channelStatus != "Normal") {
298 const KLMElectronicsChannel* electronicsChannel =
299 m_ElectronicsMap->getElectronicsChannel(channelNumber);
300 if (electronicsChannel == nullptr)
301 B2FATAL("Incomplete BKLM electronics map.");
302 if (channelStatus == "Masked") {
303 histogram->SetBinContent(i, 0);
304 if (delta != nullptr)
305 delta->SetBinContent(i, 0);
306 }
307 str = channelStatus + " channel: ";
308 // lane, axis, channel
309 str += ("L" + std::to_string(electronicsChannel->getLane()) +
310 " A" + std::to_string(electronicsChannel->getAxis()) +
311 " Ch" + std::to_string(electronicsChannel->getChannel()));
312 message_counter++;
313 if (message_counter <= m_MessageThreshold) {
314 latex.DrawLatexNDC(x, y, str.c_str());
315 y -= 0.05;
316 }
317 }
318 }
319 if (message_counter > m_MessageThreshold) {
320 std::string verbose_message = " more messages";
321 verbose_message = std::to_string(message_counter - m_MessageThreshold) + verbose_message;
322 latex.DrawLatexNDC(x, y, verbose_message.c_str());
323 }
324
325
326 // for hot/masked channels, log scale plots (reference and main)
327 if (histogram->GetMaximum()*n > histogram->Integral()*m_ThresholdForLog && average * activeModuleChannels > m_MinHitsForFlagging) {
328 histogram->SetMinimum(1);
329 canvas->SetLogy();
330 } else if (ref_histogram != nullptr) {
331 if (ref_histogram->GetMaximum()*n > ref_histogram->Integral()*m_ThresholdForLog
332 && ref_average * activeModuleChannels > m_MinHitsForFlagging) {
333 histogram->SetMinimum(1);
334 canvas->SetLogy();
335 } else
336 canvas->SetLogy(0);
337 } else
338 canvas->SetLogy(0);
339
340 canvas->Modified();
341 canvas->Update();
342
343 /* Drawing dividing lines */
344 int divisions;
345 int bin = 1;
346 double xLine;
347 // drawing lines for BKLM sectors
348 if (subdetector == 1) {
349 int shift;
350 if (index == 0) {
351 divisions = 7;
352 shift = 1;
353 } else {
356 }
357 for (int k = 0; k < divisions; k++) {
358 xLine = (histogram->GetXaxis()->GetBinLowEdge(bin) - canvas->GetX1()) / (canvas->GetX2() - canvas->GetX1());
359 m_PlaneLine.DrawLineNDC(xLine, 0.1, xLine, 0.9);
360 bin += BKLMElementNumbers::getNStrips(section, sector, k + shift, 0)
361 + BKLMElementNumbers::getNStrips(section, sector, k + shift, 1);
362 }
363 } else { // drawing lines for EKLM sectors
364 if ((section == 2) && (index == 0 || index == 1))
365 divisions = 5;
366 else
367 divisions = 4;
368 for (int k = 0; k < divisions; k++) {
369 xLine = (histogram->GetXaxis()->GetBinLowEdge(bin) - canvas->GetX1()) / (canvas->GetX2() - canvas->GetX1());
370 m_PlaneLine.DrawLineNDC(xLine, 0.1, xLine, 0.9);
372 }
373 }
374 canvas->Modified();
375 canvas->Update();
376
377}
378
380 KLMSectionNumber section, TH2F* histogram, TCanvas* canvas)
381{
382 canvas->Clear();
383 canvas->cd();
384 histogram->SetStats(false);
385 histogram->Draw("COLZ");
386 /* Draw the lines only for the backward layers. */
388 m_2DHitsLine.DrawLine(-110, 80, -110, 190);
389 m_2DHitsLine.DrawLine(-110, 190, 110, 190);
390 m_2DHitsLine.DrawLine(110, 80, 110, 190);
391 m_2DHitsLine.DrawLine(-110, 80, 110, 80);
392 }
393 canvas->Modified();
394 canvas->Update();
395}
396
398 const std::string& histName)
399{
400 TH1* histogram = findHist(m_histogramDirectoryName + "/" + histName);
401 if (histogram == nullptr) {
402 B2WARNING("KLM DQM histogram " + m_histogramDirectoryName + "/" << histName << " is not found.");
403 return;
404 }
405
406 TCanvas* canvas = findCanvas(m_histogramDirectoryName + "/c_" + histName);
407 if (canvas == nullptr) {
408 B2WARNING("KLM DQM histogram canvas " + m_histogramDirectoryName + "/c_" << histName << " is not found.");
409 return;
410 }
411
412 else {
413 canvas->Clear();
414 canvas->cd();
415 histogram->Draw();
416 /* calling on delta histogram*/
417 TH1* delta = getDelta(m_histogramDirectoryName, histName);
418 UpdateCanvas(canvas->GetName(), delta != nullptr); //keeping this for testing purposes
419 if (delta != nullptr) {
420 B2INFO("DQMHistAnalysisKLM: Time Delta Entries is " << delta->GetEntries());
421 deltaDrawer(delta, histogram, canvas);
422 }
423 }
424}
425
427 const std::string& histName)
428{
429 TH1* histogram = findHist(m_histogramDirectoryName + "/" + histName);
430 if (histogram == nullptr) {
431 B2WARNING("KLM DQM histogram " + m_histogramDirectoryName + "/" << histName << " is not found.");
432 return;
433 }
434 TCanvas* canvas = findCanvas(m_histogramDirectoryName + "/c_" + histName);
435 if (canvas == nullptr) {
436 B2WARNING("KLM DQM histogram canvas " + m_histogramDirectoryName + "/c_" << histName << " is not found.");
437 return;
438 }
439
440 histogram->Clear();
441 canvas->Clear();
442 canvas->cd();
443 if (m_MaskedChannels.size() > 0) {
444 int channelSubdetector, channelSection, channelSector;
445 int layer, plane, strip;
446 for (KLMChannelNumber channel : m_MaskedChannels) {
448 channel, &channelSubdetector, &channelSection, &channelSector,
449 &layer, &plane, &strip);
450 KLMSectorNumber sectorNumber;
451 if (channelSubdetector == KLMElementNumbers::c_BKLM)
452 sectorNumber = m_ElementNumbers->sectorNumberBKLM(channelSection, channelSector);
453 else
454 sectorNumber = m_ElementNumbers->sectorNumberEKLM(channelSection, channelSector);
455 KLMSectorNumber sectorIndex = m_SectorArrayIndex->getIndex(sectorNumber);
456 histogram->Fill(sectorIndex);
457 }
458 }
459 histogram->SetStats(false);
460 histogram->Draw();
461 canvas->Modified();
462 canvas->Update();
463}
464
466 const std::string& histName, TLatex& latex)
467{
468
469 TH1* histogram = findHist(m_histogramDirectoryName + "/" + histName);
470 if (histogram == nullptr) {
471 B2WARNING("KLM DQM histogram " + m_histogramDirectoryName + "/" << histName << " is not found.");
472 return;
473 }
474 TCanvas* canvas = findCanvas(m_histogramDirectoryName + "/c_" + histName);
475 if (canvas == nullptr) {
476 B2WARNING("KLM DQM histogram canvas " + m_histogramDirectoryName + "/c_" << histName << " is not found.");
477 return;
478 } else {
479 std::string name, alarm;
480 int moduleSubdetector, moduleSection, moduleSector, moduleLayer;
481 double xAlarm = 0.15;
482 double yAlarm = 0.8;
483 canvas->Clear();
484 canvas->cd();
485 histogram->SetStats(false);
486 histogram->Draw();
487
488 int message_counter = 0;
489 if (histName.find("bklm") != std::string::npos) {
490 /* First draw the vertical lines and the sector names. */
491 const int maximalLayer = BKLMElementNumbers::getMaximalLayerNumber();
492 for (int sector = 0; sector < BKLMElementNumbers::getMaximalSectorGlobalNumber(); ++sector) {
493 int bin = maximalLayer * sector + 1;
494 double xLine = histogram->GetXaxis()->GetBinLowEdge(bin);
495 double xText = histogram->GetXaxis()->GetBinLowEdge(bin + maximalLayer / 2);
496 double yText = gPad->GetUymin() + 0.98 * (gPad->GetUymax() - gPad->GetUymin());
497 if (sector > 0)
498 m_PlaneLine.DrawLine(xLine, gPad->GetUymin(), xLine, gPad->GetUymax());
499 name = "B";
501 name += "B";
502 else
503 name += "F";
504 name += std::to_string(sector % BKLMElementNumbers::getMaximalSectorNumber());
505 m_PlaneText.DrawText(xText, yText, name.c_str());
506 }
507 /* Then, color the canvas with red if there is a dead module
508 * and write an error message. */
509 if (m_DeadBarrelModules.size() == 0) {
512 for (KLMModuleNumber module : m_DeadBarrelModules) {
514 module, &moduleSubdetector, &moduleSection, &moduleSector, &moduleLayer);
515 alarm = "No data from " + m_ElementNumbers->getSectorDAQName(moduleSubdetector, moduleSection, moduleSector);
516 alarm += ", layer " + std::to_string(moduleLayer);
517 message_counter++;
518 if (message_counter <= m_MessageThreshold) {
519 latex.DrawLatexNDC(xAlarm, yAlarm, alarm.c_str());
520 yAlarm -= 0.05;
521 }
522 }
523 if (m_IsNullRun == false) {
525 }
526 } //end of enough statistics condition
527 else {
529 }
530 } else {
531 /* First draw the vertical lines and the sector names. */
532 const double maximalLayer = EKLMElementNumbers::getMaximalLayerGlobalNumber();
534 for (int layerGlobal = 1; layerGlobal <= maximalLayer; ++layerGlobal) {
535 int bin = maxPlane * layerGlobal + 1;
536 double xLine = histogram->GetXaxis()->GetBinLowEdge(bin);
537 double xText = histogram->GetXaxis()->GetBinLowEdge(bin - maxPlane / 2);
538 double yText = gPad->GetUymin() + 0.98 * (gPad->GetUymax() - gPad->GetUymin());
539 if (layerGlobal < maximalLayer)
540 m_PlaneLine.DrawLine(xLine, gPad->GetUymin(), xLine, gPad->GetUymax());
541 int section, layer;
543 layerGlobal, &section, &layer);
545 name = "B";
546 else
547 name = "F";
548 name += std::to_string(layer);
549 m_PlaneText.DrawText(xText, yText, name.c_str());
550 }
551 /* Then, color the canvas with red if there is a dead module
552 * and write an error message. */
553 if (m_DeadEndcapModules.size() == 0) {
556 for (KLMModuleNumber module : m_DeadEndcapModules) {
558 module, &moduleSubdetector, &moduleSection, &moduleSector, &moduleLayer);
559 alarm = "No data from " + m_ElementNumbers->getSectorDAQName(moduleSubdetector, moduleSection, moduleSector);
560 alarm += ", layer " + std::to_string(moduleLayer);
561 message_counter++;
562 if (message_counter <= m_MessageThreshold) {
563 latex.DrawLatexNDC(xAlarm, yAlarm, alarm.c_str());
564 yAlarm -= 0.05;
565 }
566 }
567 if (m_IsNullRun == false) {
569 }
570 } //end of high statistics condition
571 else {
573 }
574 }
575 if (message_counter > m_MessageThreshold) {
576 std::string verbose_string = " more messages";
577 verbose_string = std::to_string(message_counter - m_MessageThreshold) + verbose_string;
578 latex.DrawLatexNDC(xAlarm, yAlarm, verbose_string.c_str());
579 }
580 canvas->Modified();
581 canvas->Update();
582 }
583}
584
586{
587 /* If KLM is not included, stop here and return. */
588 TH1* daqInclusion = findHist(m_histogramDirectoryName + "/daq_inclusion");
589 if (not(daqInclusion == nullptr)) {
590 int isKlmIncluded = daqInclusion->GetBinContent(daqInclusion->GetXaxis()->FindBin("Yes"));
591 if (isKlmIncluded == 0)
592 return;
593 }
594 /* Make sure that the vectors are cleared at each DQM refresh. */
595 m_DeadBarrelModules.clear();
596 m_DeadEndcapModules.clear();
597 m_MaskedChannels.clear();
599 std::string str, histogramName, canvasName;
600 TLatex latex;
601 latex.SetTextColor(kRed);
602 latex.SetTextAlign(11);
604 // gathering relevant info for analyseChannelHitHistogram
605 for (KLMChannelIndex& klmSector : klmIndex) {
606 int nHistograms;
607 if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM)
608 nHistograms = 2;
609 else
610 nHistograms = 3;
611 for (int j = 0; j < nHistograms; j++) {
612 str = "strip_hits_subdetector_" +
613 std::to_string(klmSector.getSubdetector()) +
614 "_section_" + std::to_string(klmSector.getSection()) +
615 "_sector_" + std::to_string(klmSector.getSector()) +
616 "_" + std::to_string(j);
617 histogramName = m_histogramDirectoryName + "/" + str;
618 canvasName = m_histogramDirectoryName + "/c_" + str;
619 TH1* histogram = findHist(histogramName);
620 //get delta histogram (should we work with a clone instead?)
621 auto delta = getDelta("", histogramName);
622
623 if (histogram == nullptr) {
624 B2WARNING("KLM DQM histogram " << histogramName << " is not found.");
625 continue;
626 }
627 TCanvas* canvas = findCanvas(canvasName);
628 if (canvas == nullptr) {
629 B2WARNING("KLM DQM histogram canvas " << canvasName << " is not found.");
630 continue;
631 }
632 // Add this canvas that it is time to update
633 // not sure if this is interfering with the generation of some features
634 // after testing, switch condition back to delta != nullptr || histogram != nullptr
635 UpdateCanvas(canvas->GetName(), true);
637 klmSector.getSubdetector(), klmSector.getSection(),
638 klmSector.getSector(), j, histogram, delta, canvas, latex);
639
640 }
641 }
642 /* Temporary change the color palette. */
643 gStyle->SetPalette(kLightTemperature);
645 for (KLMChannelIndex& klmSection : klmIndex) {
646 KLMSubdetectorNumber subdetector = klmSection.getSubdetector();
647 if (subdetector == KLMElementNumbers::c_EKLM) {
648 KLMSubdetectorNumber section = klmSection.getSection();
649 int maximalLayerNumber = m_EklmElementNumbers->getMaximalDetectorLayerNumber(section);
650 for (int j = 1; j <= maximalLayerNumber; ++j) {
651 str = "spatial_2d_hits_subdetector_" + std::to_string(subdetector) +
652 "_section_" + std::to_string(section) +
653 "_layer_" + std::to_string(j);
654 histogramName = m_histogramDirectoryName + "/" + str;
655 canvasName = m_histogramDirectoryName + "/c_" + str;
656 TH2F* histogram = static_cast<TH2F*>(findHist(histogramName));
657 if (histogram == nullptr) {
658 B2ERROR("KLM DQM histogram " << histogramName << " is not found.");
659 continue;
660 }
661 TCanvas* canvas = findCanvas(canvasName);
662 if (canvas == nullptr) {
663 B2ERROR("KLM DQM histogram canvas " << canvasName << " is not found.");
664 continue;
665 }
666 processSpatial2DHitEndcapHistogram(section, histogram, canvas);
667 }
668 }
669 }
670 /* Reset the color palette to the default one. */
671 gStyle->SetPalette(kBird);
672 fillMaskedChannelsHistogram("masked_channels");
673 latex.SetTextColor(kBlue);
674 processPlaneHistogram("plane_bklm_phi", latex);
675 processPlaneHistogram("plane_bklm_z", latex);
676 processPlaneHistogram("plane_eklm", latex);
677
678 processTimeHistogram("time_rpc");
679 processTimeHistogram("time_scintillator_bklm");
680 processTimeHistogram("time_scintillator_eklm");
681
682 B2DEBUG(20, "Updating EPICS PVs for DQMHistAnalysisKLM");
683 // only update PVs if there's enough statistics and datasize != 0
684 // Check if it's a null run, if so, don't update EPICS PVs
685 if (m_IsNullRun) {
686 B2INFO("DQMHistAnalysisKLM: Null run detected. No PV Update.");
687 return;
688 }
689 auto* daqDataSize = findHist("DAQ/KLMDataSize");
690 double meanDAQDataSize = 0.;
691 if (daqDataSize != nullptr) {
692 meanDAQDataSize = daqDataSize->GetMean();
693 } else
694 B2WARNING("DQMHistAnalysisKLM: Cannot find KLMDataSize");
695 if ((daqDataSize != nullptr) and (meanDAQDataSize != 0.)) {
697 setEpicsPV("MaskedChannels", (double)m_MaskedChannels.size());
698 setEpicsPV("DeadBarrelModules", (double)m_DeadBarrelModules.size());
699 setEpicsPV("DeadEndcapModules", (double)m_DeadEndcapModules.size());
700 B2DEBUG(20, "DQMHistAnalysisKLM: MaskedChannels " << m_MaskedChannels.size());
701 B2DEBUG(20, "DQMHistAnalysisKLM: DeadBarrelModules " << m_DeadBarrelModules.size());
702 B2DEBUG(20, "DQMHistAnalysisKLM: DeadEndcapModules " << m_DeadEndcapModules.size());
703 }
704 } else
705 B2INFO("DQMHistAnalysisKLM: KLM Not included. No PV Update. ");
706
707
708}
static constexpr int getMaximalLayerNumber()
Get maximal layer number (1-based).
static constexpr int getMaximalSectorNumber()
Get maximal sector number (1-based).
static constexpr int getMaximalSectorGlobalNumber()
Get maximal sector global number.
static int getNStrips(int section, int sector, int layer, int plane)
Get number of strips.
TLine m_PlaneLine
TLine for boundary in plane histograms.
double m_ProcessedEvents
Number of processed events.
void initialize() override final
Initializer.
int m_MessageThreshold
Message Threshold for expert pots.
int m_ThresholdForLog
Threshold for log scale.
double m_minEntries
Minimal number of entries for delta histogram update.
int m_MinHitsForFlagging
Minimal number of hits for flagging.
const KLMElementNumbers * m_ElementNumbers
KLM element numbers.
const EKLMElementNumbers * m_EklmElementNumbers
EKLM element numbers.
void processSpatial2DHitEndcapHistogram(uint16_t section, TH2F *histogram, TCanvas *canvas)
Process spatial 2D hits histograms for endcap.
MonitoringObject * m_monObj
Monitoring object.
void deltaDrawer(TH1 *delta, TH1 *histogram, TCanvas *canvas)
Scales and draws desired delta histogram for current canvas.
void processTimeHistogram(const std::string &histName)
Process histogram containing the hit times.
std::vector< uint16_t > m_MaskedChannels
Vector of masked channels.
void terminate() override final
This method is called at the end of the event processing.
void event() override final
This method is called for each event.
int m_ThresholdForHot
Threshold for hot channels.
std::vector< uint16_t > m_DeadBarrelModules
Vector of dead barrel modules.
std::string m_histogramDirectoryName
Name of histogram directory.
double m_MinProcessedEventsForMessagesInput
Input parameter for minimal number of processed events for error messages.
TText m_PlaneText
TText for names in plane histograms.
DBObjPtr< KLMElectronicsMap > m_ElectronicsMap
Electronics map.
void endRun() override final
This method is called if the current run ends.
void processPlaneHistogram(const std::string &histName, TLatex &latex)
Process histogram containing the number of hits in plane.
int m_ThresholdForMasked
Threshold for masked channels.
const KLMSectorArrayIndex * m_SectorArrayIndex
KLM sector array index.
void beginRun() override final
Called when entering a new run.
double m_MinProcessedEventsForMessages
Minimal number of processed events for error messages.
std::vector< uint16_t > m_DeadEndcapModules
Vector of dead endcap modules.
double getProcessedEvents()
Get number of processed events.
TLine m_2DHitsLine
TLine for background region in 2d hits histograms.
const KLMChannelArrayIndex * m_ChannelArrayIndex
KLM channel array index.
bool m_IsNullRun
Run type flag for null runs.
void analyseChannelHitHistogram(int subdetector, int section, int sector, int index, TH1 *histogram, TH1 *delta, TCanvas *canvas, TLatex &latex)
Analyse channel hit histogram.
void fillMaskedChannelsHistogram(const std::string &histName)
Fill histogram containing masked channels per sector.
The base class for the histogram analysis module.
TCanvas * findCanvas(TString cname)
Find canvas by name.
static MonitoringObject * getMonitoringObject(const std::string &name)
Get MonitoringObject with given name (new object is created if non-existing)
void addDeltaPar(const std::string &dirname, const std::string &histname, HistDelta::EDeltaType t, int p, unsigned int a=1)
Add Delta histogram parameters.
void colorizeCanvas(TCanvas *canvas, EStatus status)
Helper function for Canvas colorization.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
static const std::string & getRunType(void)
Get the Run Type.
TH1 * getDelta(const std::string &fullname, int n=0, bool onlyIfUpdated=true)
Get Delta histogram.
static TH1 * findRefHist(const std::string &histname, int scaling=0, const TH1 *hist=nullptr)
Get referencehistogram from list (no other search).
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
@ c_StatusTooFew
Not enough entries/event to judge.
@ c_StatusError
Analysis result: Severe issue found.
@ c_StatusGood
Analysis result: Good.
static int getEventProcessed(void)
Get the number of processed events.
int registerEpicsPV(std::string pvname, std::string keyname="")
EPICS related Functions.
void UpdateCanvas(std::string name, bool updated=true)
Mark canvas as updated (or not)
EKLM element numbers.
static constexpr int getNStripsSector()
Get number of strips in a sector.
static constexpr int getMaximalLayerGlobalNumber()
Get maximal detector layer global number.
int getMaximalDetectorLayerNumber(int section) const
Get maximal detector layer number.
void layerNumberToElementNumbers(int layerGlobal, int *section, int *layer) const
Get element numbers by detector layer global number.
static constexpr int getMaximalSectorNumber()
Get maximal sector number.
static constexpr int getMaximalPlaneNumber()
Get maximal plane number.
KLM channel array index.
KLM channel index.
void setIndexLevel(enum IndexLevel indexLevel)
Set index level.
BKLM electronics channel.
int getChannel() const
Get channel.
uint16_t getNumber(uint16_t index) const
Get element number.
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 channelNumber(int subdetector, int section, int sector, int layer, int plane, int strip) const
Get channel number.
void channelNumberToElementNumbers(KLMChannelNumber channel, int *subdetector, int *section, int *sector, int *layer, int *plane, int *strip) const
Get element numbers by channel number.
unsigned int getNChannelsModule(KLMModuleNumber module) const
Get number of channels in module.
void moduleNumberToElementNumbers(KLMModuleNumber module, int *subdetector, int *section, int *sector, int *layer) const
Get element numbers by module number.
KLMModuleNumber moduleNumber(int subdetector, int section, int sector, int layer) const
Get module number.
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.
KLM sector array index.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setVariable(const std::string &var, float val, float upErr=-1., float dwErr=-1)
set value to float variable (new variable is made if not yet existing)
Class to store variables with their name which were sent to the logging service.
TH2 * moduleHitMap(TH1 *hitMap, int moduleID)
Make hit map in HAPD view (12*12 channels).
Definition: hitMapMaker.cc:38
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 KLMModuleNumber
Module number.
uint16_t KLMSectionNumber
Section number.
Abstract base class for different kinds of events.