Belle II Software  release-08-01-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  // read thresholds from EPICS
57  requestLimitsFromEpicsPVs("abortRate", dummy_lowerAlarm, dummy_lowerWarn, dummy_upperWarn, buffThreshold);
58  requestLimitsFromEpicsPVs("minNoEvents", dummy_lowerAlarm, buffMinEvents, dummy_upperWarn, dummy_upperAlarm);
59 
60  if (!std::isnan(buffThreshold)) {
61  B2INFO(getName() << ": Setting failure rate threshold from EPICS. New failureRateThreshold " << buffThreshold);
62  m_failureRateThreshold = buffThreshold;
63  }
64  if (!std::isnan(buffMinEvents)) {
65  B2INFO(getName() << ": Setting min number of events threshold from EPICS. New minNoEvents " << buffMinEvents);
66  m_statThreshold = buffMinEvents;
67  }
68 
69 }
70 
71 
73 {
74 
75  gROOT->cd();
76  m_cAbortRate = new TCanvas("TrackingAnalysis/c_AbortRate");
77  m_cAbortRateHER = new TCanvas("TrackingAnalysis/c_AbortRateHER");
78  m_cAbortRateLER = new TCanvas("TrackingAnalysis/c_AbortRateLER");
79 
80  // add MonitoringObject
81  m_monObj = getMonitoringObject("trackingHLT");
82 
83  // register the PVs for setting thresholds
84  registerEpicsPV("TRACKING:minNoEvents", "minNoEvents", false);
85 
86  // variables to be monitored via EPICS
87  registerEpicsPV("trackingHLT:nTracksPerEvent", "nTracksPerEvent", false);
88  registerEpicsPV("trackingHLT:nVXDTracksPerEvent", "nVXDTracksPerEvent", false);
89  registerEpicsPV("trackingHLT:nCDCTracksPerEvent", "nCDCTracksPerEvent", false);
90  registerEpicsPV("trackingHLT:nVXDCDCTracksPerEvent", "nVXDCDCTracksPerEvent", false);
91  registerEpicsPV("trackingHLT:abortRate", "abortRate", false);
92 
93 }
94 
96 {
97 
98  //check Tracking Abort Rate
99  TH1* hAbort = findHist("TrackingHLTDQM/NumberTrackingErrorFlags");
100  if (hAbort != nullptr) {
101 
102  bool hasError = false;
103  int nEvents = hAbort->GetEntries();
104  double abortRate = hAbort->GetMean();
105  hAbort->SetTitle(Form("Fraction of Events in which Tracking aborts = %.4f %%", abortRate * 100));
106 
107  if (nEvents >= m_statThreshold) {
108  m_monObj->setVariable("abortRate", abortRate);
109  setEpicsPV("abortRate", abortRate);
110  }
111 
112  //check if number of errors is above the allowed limit
113  if (abortRate > m_failureRateThreshold)
114  hasError = true;
115 
116  m_cAbortRate->cd();
117  hAbort->Draw();
118 
119  if (nEvents < m_statThreshold) colorizeCanvas(m_cAbortRate, EStatus::c_StatusTooFew);
120  else if (hasError) colorizeCanvas(m_cAbortRate, EStatus::c_StatusError);
121  else colorizeCanvas(m_cAbortRate, EStatus::c_StatusGood);
122 
123  } else { // histogram not found
124  B2WARNING("Histogram TrackingHLTDQM/NumberTrackingErrorFlags from Tracking DQM not found!");
125  m_cAbortRate->SetFillColor(kGray);
126  }
127 
128  m_cAbortRate->Modified();
129  m_cAbortRate->Update();
130 
131 
132  if (m_printCanvas)
133  m_cAbortRate->Print("c_AbortRate.pdf");
134 
135 
136  // check tracking abort rate VS time after last HER injection and time within a beam cycle HER
137  TH2F* hAbortHER = dynamic_cast<TH2F*>(findHist("TrackingHLTDQM/TrkAbortVsTimeHER"));
138  TH2F* hAllHER = dynamic_cast<TH2F*>(findHist("TrackingHLTDQM/allEvtsVsTimeHER"));
139  if (hAbortHER != nullptr && hAllHER != nullptr) {
140 
141  TH2F* hAbortRateHER = new TH2F(*hAbortHER);
142 
143  for (int i = 0; i < hAbortRateHER->GetXaxis()->GetNbins(); i++)
144  for (int j = 0; j < hAbortRateHER->GetYaxis()->GetNbins(); j++) {
145  int den = hAllHER->GetBinContent(i + 1, j + 1);
146  int num = hAbortHER->GetBinContent(i + 1, j + 1);
147 
148  if (den > 0) hAbortRateHER->SetBinContent(i + 1, j + 1, num * 1. / den);
149  else hAbortRateHER->SetBinContent(i + 1, j + 1, 0);
150  }
151 
152  m_cAbortRateHER->cd();
153  m_cAbortRateHER->SetFillColor(kWhite);
154  hAbortRateHER->SetTitle("Fraction of Events with Tracking Aborts vs HER injection");
155  hAbortRateHER->GetZaxis()->SetTitle("Fraction of events / bin");
156  hAbortRateHER->Draw("colz");
157 
158 
159  } else { // histogram not found
160  B2WARNING("Histograms TrackingHLTDQM/TrkAbortVsTimeHER or allEvtsVsTimeHER from Tracking DQM not found!");
161  m_cAbortRateHER->SetFillColor(kGray);
162  }
163  if (m_printCanvas)
164  m_cAbortRateHER->Print("c_AbortRateHER.pdf");
165 
166 
167  // check tracking abort rate VS time after last LER injection and time within a beam cycle LER
168  TH2F* hAbortLER = dynamic_cast<TH2F*>(findHist("TrackingHLTDQM/TrkAbortVsTimeLER"));
169  TH2F* hAllLER = dynamic_cast<TH2F*>(findHist("TrackingHLTDQM/allEvtsVsTimeLER"));
170  if (hAbortLER != nullptr && hAllLER != nullptr) {
171 
172  TH2F* hAbortRateLER = new TH2F(*hAbortLER);
173 
174  for (int i = 0; i < hAbortRateLER->GetXaxis()->GetNbins(); i++)
175  for (int j = 0; j < hAbortRateLER->GetYaxis()->GetNbins(); j++) {
176  int den = hAllLER->GetBinContent(i + 1, j + 1);
177  int num = hAbortLER->GetBinContent(i + 1, j + 1);
178 
179  if (den > 0) hAbortRateLER->SetBinContent(i + 1, j + 1, num * 1. / den);
180  else hAbortRateLER->SetBinContent(i + 1, j + 1, 0);
181  }
182 
183  m_cAbortRateLER->cd();
184  m_cAbortRateLER->SetFillColor(kWhite);
185  hAbortRateLER->SetTitle("Fraction of Events with Tracking Aborts vs LER injection");
186  hAbortRateLER->GetZaxis()->SetTitle("Fraction of events / bin");
187  hAbortRateLER->Draw("colz");
188 
189  } else { // histogram not found
190  B2WARNING("Histograms TrackingHLTDQM/TrkAbortVsTimeLER or allEvtsVsTimeLER from Tracking DQM not found!");
191  m_cAbortRateLER->SetFillColor(kGray);
192  }
193  if (m_printCanvas)
194  m_cAbortRateLER->Print("c_AbortRateLER.pdf");
195 
196  // add average number of tracks per event to Mirabelle
197  TH1* hnTracks = findHist("TrackingHLTDQM/NoOfTracks");
198  if (hnTracks != nullptr && hnTracks->GetEntries() >= m_statThreshold) {
199  double averageNTracks = hnTracks->GetMean();
200  m_monObj->setVariable("nTracksPerEvent", averageNTracks);
201  setEpicsPV("nTracksPerEvent", averageNTracks);
202  }
203 
204  TH1* hnVXDTracks = findHist("TrackingHLTDQM/NoOfTracksInVXDOnly");
205  if (hnVXDTracks != nullptr && hnVXDTracks->GetEntries() >= m_statThreshold) {
206  double averageNVXDTracks = hnVXDTracks->GetMean();
207  m_monObj->setVariable("nVXDTracksPerEvent", averageNVXDTracks);
208  setEpicsPV("nVXDTracksPerEvent", averageNVXDTracks);
209  }
210 
211  TH1* hnCDCTracks = findHist("TrackingHLTDQM/NoOfTracksInCDCOnly");
212  if (hnCDCTracks != nullptr && hnCDCTracks->GetEntries() >= m_statThreshold) {
213  double averageNCDCTracks = hnCDCTracks->GetMean();
214  m_monObj->setVariable("nCDCTracksPerEvent", averageNCDCTracks);
215  setEpicsPV("nCDCTracksPerEvent", averageNCDCTracks);
216  }
217 
218  TH1* hnVXDCDCTracks = findHist("TrackingHLTDQM/NoOfTracksInVXDCDC");
219  if (hnVXDCDCTracks != nullptr && hnVXDCDCTracks->GetEntries() >= m_statThreshold) {
220  double averageNVXDCDCTracks = hnVXDCDCTracks->GetMean();
221  m_monObj->setVariable("nVXDCDCTracksPerEvent", averageNVXDCDCTracks);
222  setEpicsPV("nVXDCDCTracksPerEvent", averageNVXDCDCTracks);
223  }
224 
225 }
226 
227 
228 
The base class for the histogram analysis module.
int registerEpicsPV(std::string pvname, std::string keyname="", bool update_pvs=true)
EPICS related Functions.
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).
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
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.