Belle II Software  release-08-00-10
DQMHistAnalysisTrackingHLT.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 #include <dqm/analysis/modules/DQMHistAnalysisTrackingHLT.h>
10 
11 #include <TROOT.h>
12 #include <TString.h>
13 #include <TH1F.h>
14 #include <TH2F.h>
15 
16 #include <TMath.h>
17 #include <iostream>
18 
19 using namespace std;
20 using namespace Belle2;
21 
22 //-----------------------------------------------------------------
23 // Register the Module
24 //-----------------------------------------------------------------
25 REG_MODULE(DQMHistAnalysisTrackingHLT);
26 
27 //-----------------------------------------------------------------
28 // Implementation
29 //-----------------------------------------------------------------
30 
31 DQMHistAnalysisTrackingHLTModule::DQMHistAnalysisTrackingHLTModule()
33 {
34 
35  setDescription("DQM Analysis Module of the Tracking HLT Plots.");
36 
37  addParam("failureRateThreshold", m_failureRateThreshold,
38  "Maximum Fraction of Events in which Tracking Aborts before turning Canvas to Red. Will be taken from Epics by default, \
39  this value is only taken if Epics is not available!", double(m_failureRateThreshold));
40  addParam("minNoEvents", m_statThreshold,
41  "Minimum Number of Events before scaring CR shifters. Will be taken from Epics by default, \
42  this value is only taken if Epics is not available!", int(m_statThreshold));
43  addParam("printCanvas", m_printCanvas, "if True prints pdf of the analysis canvas", bool(m_printCanvas));
44 
45 }
46 
47 
48 
50 {
51  // get the abort rate and statThreshold from epics
52  double buffThreshold(NAN);
53  double buffMinEvents(NAN);
54  double dummy_lowerAlarm, dummy_lowerWarn, dummy_upperWarn, dummy_upperAlarm;
55 
56  requestLimitsFromEpicsPVs("failureRateThreshold", dummy_lowerAlarm, dummy_lowerWarn, buffThreshold, dummy_upperAlarm);
57  requestLimitsFromEpicsPVs("minNoEvents", dummy_lowerAlarm, buffMinEvents, dummy_upperWarn, dummy_upperAlarm);
58 
59  if (!std::isnan(buffThreshold)) {
60  B2INFO(getName() << ": Setting failure rate threshold from EPICS. New failureRateThreshold " << buffThreshold);
61  m_failureRateThreshold = buffThreshold;
62  }
63  if (!std::isnan(buffMinEvents)) {
64  B2INFO(getName() << ": Setting min number of events threshold from EPICS. New minNoEvents " << buffMinEvents);
65  m_statThreshold = buffMinEvents;
66  }
67 }
68 
69 
71 {
72 
73  gROOT->cd();
74  m_cAbortRate = new TCanvas("TrackingAnalysis/c_AbortRate");
75  m_cAbortRateHER = new TCanvas("TrackingAnalysis/c_AbortRateHER");
76  m_cAbortRateLER = new TCanvas("TrackingAnalysis/c_AbortRateLER");
77 
78  // add MonitoringObject
79  m_monObj = getMonitoringObject("trackingHLT");
80 
81  // register the PVs for setting thresholds
82  registerEpicsPV("TRACKING:failureRateThreshold", "failureRateThreshold");
83  registerEpicsPV("TRACKING:minNoEvents", "minNoEvents");
84 }
85 
87 {
88 
89  //check Tracking Abort Rate
90  TH1* hAbort = findHist("TrackingHLTDQM/NumberTrackingErrorFlags");
91  if (hAbort != nullptr) {
92 
93  bool hasError = false;
94  int nEvents = hAbort->GetEntries();
95  double abortRate = hAbort->GetMean();
96  hAbort->SetTitle(Form("Fraction of Events in which Tracking aborts = %.4f %%", abortRate * 100));
97 
98  m_monObj->setVariable("abortRate", abortRate);
99  //check if number of errors is above the allowed limit
100  if (abortRate > m_failureRateThreshold)
101  hasError = true;
102 
103  if (nEvents < m_statThreshold) m_cAbortRate->SetFillColor(kGray);
104  else if (hasError) m_cAbortRate->SetFillColor(kRed);
105  else m_cAbortRate->SetFillColor(kGreen);
106  m_cAbortRate->SetFrameFillColor(10);
107 
108  m_cAbortRate->cd();
109  hAbort->Draw();
110 
111  } else { // histogram not found
112  B2WARNING("Histogram TrackingHLTDQM/NumberTrackingErrorFlags from Tracking DQM not found!");
113  m_cAbortRate->SetFillColor(kGray);
114  }
115 
116  m_cAbortRate->Modified();
117  m_cAbortRate->Update();
118 
119 
120  if (m_printCanvas)
121  m_cAbortRate->Print("c_AbortRate.pdf");
122 
123 
124  // check tracking abort rate VS time after last HER injection and time within a beam cycle HER
125  TH2F* hAbortHER = dynamic_cast<TH2F*>(findHist("TrackingHLTDQM/TrkAbortVsTimeHER"));
126  TH2F* hAllHER = dynamic_cast<TH2F*>(findHist("TrackingHLTDQM/allEvtsVsTimeHER"));
127  if (hAbortHER != nullptr && hAllHER != nullptr) {
128 
129  TH2F* hAbortRateHER = new TH2F(*hAbortHER);
130 
131  for (int i = 0; i < hAbortRateHER->GetXaxis()->GetNbins(); i++)
132  for (int j = 0; j < hAbortRateHER->GetYaxis()->GetNbins(); j++) {
133  int den = hAllHER->GetBinContent(i + 1, j + 1);
134  int num = hAbortHER->GetBinContent(i + 1, j + 1);
135 
136  if (den > 0) hAbortRateHER->SetBinContent(i + 1, j + 1, num * 1. / den);
137  else hAbortRateHER->SetBinContent(i + 1, j + 1, 0);
138  }
139 
140  m_cAbortRateHER->cd();
141  hAbortRateHER->Draw("colz");
142 
143 
144  } else { // histogram not found
145  B2WARNING("Histograms TrackingHLTDQM/TrkAbortVsTimeHER or allEvtsVsTimeHER from Tracking DQM not found!");
146  m_cAbortRateHER->SetFillColor(kGray);
147  }
148  if (m_printCanvas)
149  m_cAbortRateHER->Print("c_AbortRateHER.pdf");
150 
151 
152  // check tracking abort rate VS time after last LER injection and time within a beam cycle LER
153  TH2F* hAbortLER = dynamic_cast<TH2F*>(findHist("TrackingHLTDQM/TrkAbortVsTimeLER"));
154  TH2F* hAllLER = dynamic_cast<TH2F*>(findHist("TrackingHLTDQM/allEvtsVsTimeLER"));
155  if (hAbortLER != nullptr && hAllLER != nullptr) {
156 
157  TH2F* hAbortRateLER = new TH2F(*hAbortLER);
158 
159  for (int i = 0; i < hAbortRateLER->GetXaxis()->GetNbins(); i++)
160  for (int j = 0; j < hAbortRateLER->GetYaxis()->GetNbins(); j++) {
161  int den = hAllLER->GetBinContent(i + 1, j + 1);
162  int num = hAbortLER->GetBinContent(i + 1, j + 1);
163 
164  if (den > 0) hAbortRateLER->SetBinContent(i + 1, j + 1, num * 1. / den);
165  else hAbortRateLER->SetBinContent(i + 1, j + 1, 0);
166  }
167 
168  m_cAbortRateLER->cd();
169  hAbortRateLER->Draw("colz");
170 
171  } else { // histogram not found
172  B2WARNING("Histograms TrackingHLTDQM/TrkAbortVsTimeLER or allEvtsVsTimeLER from Tracking DQM not found!");
173  m_cAbortRateLER->SetFillColor(kGray);
174  }
175  if (m_printCanvas)
176  m_cAbortRateLER->Print("c_AbortRateLER.pdf");
177 
178  // add average number of tracks per event to Mirabelle
179  TH1* hnTracks = findHist("TrackingHLTDQM/NoOfTracks");
180  if (hnTracks != nullptr) {
181  double averageNTracks = hnTracks->GetMean();
182  m_monObj->setVariable("nTracksPerEvent", averageNTracks);
183  }
184 
185  TH1* hnVXDTracks = findHist("TrackingHLTDQM/NoOfTracksInVXDOnly");
186  if (hnVXDTracks != nullptr) {
187  double averageNVXDTracks = hnVXDTracks->GetMean();
188  m_monObj->setVariable("nVXDTracksPerEvent", averageNVXDTracks);
189  }
190 
191  TH1* hnCDCTracks = findHist("TrackingHLTDQM/NoOfTracksInCDCOnly");
192  if (hnCDCTracks != nullptr) {
193  double averageNCDCTracks = hnCDCTracks->GetMean();
194  m_monObj->setVariable("nCDCTracksPerEvent", averageNCDCTracks);
195  }
196 
197  TH1* hnVXDCDCTracks = findHist("TrackingHLTDQM/NoOfTracksInVXDCDC");
198  if (hnVXDCDCTracks != nullptr) {
199  double averageNVXDCDCTracks = hnVXDCDCTracks->GetMean();
200  m_monObj->setVariable("nVXDCDCTracksPerEvent", averageNVXDCDCTracks);
201  }
202 
203 }
204 
205 
206 
The base class for the histogram analysis module.
int registerEpicsPV(std::string pvname, std::string keyname="", bool update_pvs=true)
EPICS related Functions.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
static MonitoringObject * getMonitoringObject(const std::string &histname)
Get MonitoringObject with given name (new object is created if non-existing)
bool requestLimitsFromEpicsPVs(chid id, double &lowerAlarm, double &lowerWarn, double &upperWarn, double &upperAlarm)
Get Alarm Limits from EPICS PV.
void initialize() override final
Module function initialize.
int m_statThreshold
minimal number of events to judge
double m_failureRateThreshold
above this rate, there is maybe a problem?
MonitoringObject * m_monObj
Monitoring Object to be produced by this module, which contain defined canvases and monitoring variab...
TCanvas * m_cAbortRate
canvas for the abort rate plot
void event() override final
Module function event.
TCanvas * m_cAbortRateHER
canvas for the 2D abort rate plot for HER
TCanvas * m_cAbortRateLER
canvas for the 2D abort rate plot for LER
bool m_printCanvas
if true print the pdf of the canvases
void beginRun() override final
Module function doing stuff at beginning of a run.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
const std::string & getName() const
Returns the name of the module.
Definition: Module.h:187
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)
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
Abstract base class for different kinds of events.