Belle II Software  release-06-01-15
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 
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", double(m_failureRateThreshold));
39  addParam("minNoEvents", m_statThreshold, "Minimum Number of Events before scaring CR shifters", int(m_statThreshold));
40  addParam("printCanvas", m_printCanvas, "if True prints pdf of the analysis canvas", bool(m_printCanvas));
41 
42 }
43 
44 void DQMHistAnalysisTrackingHLTModule::initialize()
45 {
46 
47  gROOT->cd();
48  m_cAbortRate = new TCanvas("TrackingAnalysis/c_AbortRate");
49  m_cAbortRateHER = new TCanvas("TrackingAnalysis/c_AbortRateHER");
50  m_cAbortRateLER = new TCanvas("TrackingAnalysis/c_AbortRateLER");
51 
52  // add MonitoringObject
53  m_monObj = getMonitoringObject("trackingHLT");
54 
55 }
56 
57 void DQMHistAnalysisTrackingHLTModule::event()
58 {
59 
60  //check Tracking Abort Rate
61 
62  TH1* hAbort = findHist("TrackingHLTDQM/NumberTrackingErrorFlags");
63  if (hAbort != nullptr) {
64 
65  bool hasError = false;
66  int nEvents = hAbort->GetEntries();
67  double abortRate = hAbort->GetMean();
68  hAbort->SetTitle(Form("Fraction of Events in which Tracking aborts = %.4f %%", abortRate * 100));
69  m_monObj->setVariable("abortRate", abortRate);
70  //check if number of errors is above the allowed limit
71  if (abortRate > m_failureRateThreshold)
72  hasError = true;
73 
74  if (nEvents < m_statThreshold) m_cAbortRate->SetFillColor(kGray);
75  else if (hasError) m_cAbortRate->SetFillColor(kRed);
76  else m_cAbortRate->SetFillColor(kGreen);
77  m_cAbortRate->SetFrameFillColor(10);
78 
79  m_cAbortRate->cd();
80  hAbort->Draw();
81 
82  } else { // histogram not found
83  B2WARNING("Histogram TrackingHLTDQM/NumberTrackingErrorFlags from Tracking DQM not found!");
84  m_cAbortRate->SetFillColor(kGray);
85  }
86 
87  m_cAbortRate->Modified();
88  m_cAbortRate->Update();
89 
90 
91  if (m_printCanvas)
92  m_cAbortRate->Print("c_AbortRate.pdf");
93 
94 
95  // check tracking abort rate VS time after last HER injection and time within a beam cycle HER
96  TH2F* hAbortHER = dynamic_cast<TH2F*>(findHist("TrackingHLTDQM/TrkAbortVsTimeHER"));
97  TH2F* hAllHER = dynamic_cast<TH2F*>(findHist("TrackingHLTDQM/allEvtsVsTimeHER"));
98  if (hAbortHER != nullptr && hAllHER != nullptr) {
99 
100  TH2F* hAbortRateHER = new TH2F(*hAbortHER);
101 
102  for (int i = 0; i < hAbortRateHER->GetXaxis()->GetNbins(); i++)
103  for (int j = 0; j < hAbortRateHER->GetYaxis()->GetNbins(); j++) {
104  int den = hAllHER->GetBinContent(i + 1, j + 1);
105  int num = hAbortHER->GetBinContent(i + 1, j + 1);
106 
107  if (den > 0) hAbortRateHER->SetBinContent(i + 1, j + 1, num * 1. / den);
108  else hAbortRateHER->SetBinContent(i + 1, j + 1, 0);
109  }
110 
111  m_cAbortRateHER->cd();
112  hAbortRateHER->Draw("colz");
113 
114 
115  } else { // histogram not found
116  B2WARNING("Histograms TrackingHLTDQM/TrkAbortVsTimeHER or allEvtsVsTimeHER from Tracking DQM not found!");
117  m_cAbortRateHER->SetFillColor(kGray);
118  }
119  if (m_printCanvas)
120  m_cAbortRateHER->Print("c_AbortRateHER.pdf");
121 
122 
123  // check tracking abort rate VS time after last LER injection and time within a beam cycle LER
124  TH2F* hAbortLER = dynamic_cast<TH2F*>(findHist("TrackingHLTDQM/TrkAbortVsTimeLER"));
125  TH2F* hAllLER = dynamic_cast<TH2F*>(findHist("TrackingHLTDQM/allEvtsVsTimeLER"));
126  if (hAbortLER != nullptr && hAllLER != nullptr) {
127 
128  TH2F* hAbortRateLER = new TH2F(*hAbortLER);
129 
130  for (int i = 0; i < hAbortRateLER->GetXaxis()->GetNbins(); i++)
131  for (int j = 0; j < hAbortRateLER->GetYaxis()->GetNbins(); j++) {
132  int den = hAllLER->GetBinContent(i + 1, j + 1);
133  int num = hAbortLER->GetBinContent(i + 1, j + 1);
134 
135  if (den > 0) hAbortRateLER->SetBinContent(i + 1, j + 1, num * 1. / den);
136  else hAbortRateLER->SetBinContent(i + 1, j + 1, 0);
137  }
138 
139  m_cAbortRateLER->cd();
140  hAbortRateLER->Draw("colz");
141 
142  } else { // histogram not found
143  B2WARNING("Histograms TrackingHLTDQM/TrkAbortVsTimeLER or allEvtsVsTimeLER from Tracking DQM not found!");
144  m_cAbortRateLER->SetFillColor(kGray);
145  }
146  if (m_printCanvas)
147  m_cAbortRateLER->Print("c_AbortRateLER.pdf");
148 
149  // add average number of tracks per event to Mirabelle
150  TH1* hnTracks = findHist("TrackingHLTDQM/NoOfTracks");
151  if (hnTracks != nullptr) {
152  double averageNTracks = hnTracks->GetMean();
153  m_monObj->setVariable("nTracksPerEvent", averageNTracks);
154  }
155 
156  TH1* hnVXDTracks = findHist("TrackingHLTDQM/NoOfTracksInVXDOnly");
157  if (hnVXDTracks != nullptr) {
158  double averageNVXDTracks = hnVXDTracks->GetMean();
159  m_monObj->setVariable("nVXDTracksPerEvent", averageNVXDTracks);
160  }
161 
162  TH1* hnCDCTracks = findHist("TrackingHLTDQM/NoOfTracksInCDCOnly");
163  if (hnCDCTracks != nullptr) {
164  double averageNCDCTracks = hnCDCTracks->GetMean();
165  m_monObj->setVariable("nCDCTracksPerEvent", averageNCDCTracks);
166  }
167 
168  TH1* hnVXDCDCTracks = findHist("TrackingHLTDQM/NoOfTracksInVXDCDC");
169  if (hnVXDCDCTracks != nullptr) {
170  double averageNVXDCDCTracks = hnVXDCDCTracks->GetMean();
171  m_monObj->setVariable("nVXDCDCTracksPerEvent", averageNVXDCDCTracks);
172  }
173 
174 }
The base class for the histogram analysis module.
#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.