ATLAS Offline Software
Loading...
Searching...
No Matches
MonitoringFile_HLTMETDQFlagSummary.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3 */
4
5/* HLTMET Post Processing Method: In the HLT/MET[Mon/Mon_allCells/Mon_FEB]/DQPlots directory, there are
6 * three summary histograms: L2_MET-status, EF_MET_status and and compN_compSumEt_lin
7 * and lumi-block based histogram: trmet_lbn_flag
8 * This post-processing method fills these histograms appropriately
9 * Author : Venkatesh Kaushik <venkat.kaushik@cern.ch>
10 * Date : Mar 2010
11 */
12
14
15#include <iostream>
16#include <iomanip>
17#include <algorithm>
18#include <fstream>
19#include <cmath>
20#include <cstdlib>
21#include <sstream>
22#include <vector>
23#include <utility>
24#include <map>
25#include <string>
26
27#include "TH1F.h"
28#include "TH1I.h"
29#include "TH2F.h"
30#include "TFile.h"
31#include "TClass.h"
32#include "TKey.h"
33#include "TMath.h"
34#include "TF1.h"
35#include "TTree.h"
36#include "TBranch.h"
37
38namespace dqutils {
39 //---------------------------------------------------------------------------------------------------
40 // Method to obtain summary DQFlags
41 //---------------------------------------------------------------------------------------------------
42 void MonitoringFile::HLTMETDQFlagSummary(TFile* f, std::string& run_dir) {
43 //bool dbgLevel = false;
44
45 //if (dbgLevel) std::cout << "--> HLTMETDQFlagSummary: Updating histograms in HLT/METMon*/DQPlots " << std::endl;
46
47 f->cd("/");
48 TIter next_run(f->GetListOfKeys());
49 TKey* key_run(0);
50 while ((key_run = dynamic_cast<TKey*> (next_run())) != 0) {
51 if (!key_run->IsFolder()) continue;
52 run_dir = key_run->GetName();
53 if (run_dir.find("run") == std::string::npos) {
54 continue;
55 }
56
57 TDirectory* tdir_run = dynamic_cast<TDirectory*>(key_run);
58
59
60
61 // begin HLTMET
62 // note 1: prefix all dirs and hists with '/'
63 // note 2: missing dir => return
64 // note 3: missing hist => continue
65 static const std::string hlt_top = run_dir + "/HLT"; // toplevel
66 static const std::string met_efdir = "/EFMissingET_Fex"; // EF dir
67 static const std::string met_l2dir = "/L2MissingET_Fex"; // L2 dir
68 static const std::string dqflag_dir = "/DQFlags"; // DQ flags dir
69
70 std::vector<std::string> met_fexs, met_l2hists, met_efhists;
71 // expect the following fex dirs
72 met_fexs.push_back("/METMon");
73 //met_fexs.push_back("/METMon_FEB");
74 //met_fexs.push_back("/METMon_allCells");
75
76 // expect the following histograms
77 met_l2hists.push_back("/L2_MET_status");
78 met_l2hists.push_back("/L2_MEx_log");
79 met_l2hists.push_back("/L2_MEy_log");
80 met_l2hists.push_back("/L2_MET_log");
81 met_l2hists.push_back("/L2_SumEt_log");
82
83 met_efhists.push_back("/EF_MET_status");
84 met_efhists.push_back("/EF_MET_lin");
85 met_efhists.push_back("/EF_MET_lin1");
86 met_efhists.push_back("/EF_SumEt_lin");
87 met_efhists.push_back("/EF_MET_phi");
88 met_efhists.push_back("/EF_MET_phi1");
89 met_efhists.push_back("/EF_MEx_log");
90 met_efhists.push_back("/EF_MEy_log");
91 met_efhists.push_back("/compN_compSumEt_lin");
92 met_efhists.push_back("/compN_EF_MET_status");
93
94 std::vector<std::string> lbnDirs;
95 static const std::string lbn_dqhist = "/trmet_lbn_flag";
96 size_t lbn_range = HLTMETGetDQLBNRange(tdir_run, lbnDirs);
97
98 //std::cout << "lbn_range = " << lbn_range << std::endl;
99
100 // Get EF/L2 status histograms for default Fex
101 for (std::vector<std::string>::iterator itFex = met_fexs.begin(); itFex != met_fexs.end(); ++itFex) {
102 std::string theL2Path = hlt_top + (*itFex) + met_l2dir;
103 std::string theEFPath = hlt_top + (*itFex) + met_efdir;
104 std::string theDQPath = hlt_top + (*itFex) + dqflag_dir;
105
106 TDirectory* dirl2 = f->GetDirectory(theL2Path.c_str());
107 if (!dirl2) {
108 std::cerr << "--> HLTMETDQFlagSummary: directory " << theL2Path << " not found" << std::endl;
109 return;
110 }
111 TDirectory* diref = f->GetDirectory(theEFPath.c_str());
112 if (!diref) {
113 std::cerr << "--> HLTMETDQFlagSummary: directory " << theEFPath << " not found" << std::endl;
114 return;
115 }
116 TDirectory* dirdq = f->GetDirectory(theDQPath.c_str());
117 if (!dirdq) {
118 std::cerr << "--> HLTMETDQFlagSummary: directory " << theDQPath << " not found" << std::endl;
119 return;
120 }
121 // loop over L2 hists and copy to DQFlags
122 for (std::vector<std::string>::iterator itHist = met_l2hists.begin(); itHist != met_l2hists.end(); ++itHist) {
123 std::string histL2 = (theL2Path + *itHist);
124 std::string histL2C = (theDQPath + *itHist);
125 TH1* hl2(0), *hl2c(0);
126 if (!f->Get(histL2.c_str())) {
127 //if (dbgLevel) std::cerr << "--> HLTMETDQFlagSummary: no histogram " << histL2 << std::endl;
128 continue;
129 }
130 if (!f->Get(histL2C.c_str())) {
131 //if (dbgLevel) std::cerr << "--> HLTMETDQFlagSummary: no histogram " << histL2C << std::endl;
132 continue;
133 }
134 hl2 = (TH1*) (f->Get(histL2.c_str()));
135 // check if MET_status histogram is emtpy
136 bool save_status_hist = false;
137 if (*itHist == "/L2_MET_status") {
138 if (hl2->GetMaximum() < 1.e-3) {
139 save_status_hist = true;
140 for (int ibin = 1; ibin <= hl2->GetNbinsX(); ibin++) {
141 hl2->SetBinContent(ibin, 0.05);
142 }
143 }
144 }
145 hl2c = (TH1*) (f->Get(histL2C.c_str()));
146 hl2c->Reset();
147 hl2c->Add(hl2);
148 //hl2c->Sumw2();
149 //std::cout << "--> HLTMETDQFlagSummary: added histograms " << histL2 << "\t and " << histL2C << std::endl;
150 dirdq->cd();
151 hl2c->Write("", TObject::kOverwrite);
152 if (save_status_hist) {
153 if (f->cd(theL2Path.c_str())) hl2->Write("", TObject::kOverwrite);
154 }
155 } // end loop over l2 hists
156
157 // loop over EF hists and copy to DQFlags
158 for (std::vector<std::string>::iterator itHist = met_efhists.begin(); itHist != met_efhists.end(); ++itHist) {
159 std::string histEF = (theEFPath + *itHist);
160 std::string histEFC = (theDQPath + *itHist);
161 TH1* hef(0), *hefc(0);
162 if (!f->Get(histEF.c_str())) {
163 //if (dbgLevel) std::cerr << "--> HLTMETDQFlagSummary: no histogram " << histEF << std::endl;
164 continue;
165 }
166 if (!f->Get(histEFC.c_str())) {
167 //if (dbgLevel) std::cerr << "--> HLTMETDQFlagSummary: no histogram " << histEFC << std::endl;
168 continue;
169 }
170 hef = (TH1*) (f->Get(histEF.c_str()));
171 bool save_status_hist = false;
172 if (*itHist == "/EF_MET_status") {
173 if (hef->GetMaximum() < 1.e-3) {
174 save_status_hist = true;
175 for (int ibin = 1; ibin <= hef->GetNbinsX(); ibin++) {
176 hef->SetBinContent(ibin, 0.05);
177 }
178 }
179 }
180 hefc = (TH1*) (f->Get(histEFC.c_str()));
181 hefc->Reset();
182 hefc->Add(hef);
183 //hefc->Sumw2();
184 dirdq->cd();
185 if (save_status_hist) {
186 if (f->cd(theEFPath.c_str())) hef->Write("", TObject::kOverwrite);
187 }
188 hefc->Write("", TObject::kOverwrite);
189 } // end loop over ef hists
190
191 // resize lbn based histogram
192 if (lbn_range > 0) {
193 std::string histLBN = theDQPath + lbn_dqhist;
194 TH1* hlb(0);
195 if (!f->Get(histLBN.c_str())) {
196 //if (dbgLevel) std::cerr << "--> HLTMETDQFlagSummary: no histogram " << histLBN << std::endl;
197 continue;
198 }
199 hlb = (TH1*) (f->Get(histLBN.c_str()));
200 unsigned int nbinx = (unsigned int) hlb->GetNbinsX(), nbiny = (unsigned int) hlb->GetNbinsY(), k = 1;
201 if (nbinx != nbiny) continue;
202
203 if (lbn_range < nbinx) {
204 hlb->GetXaxis()->SetRangeUser(0, Double_t(lbn_range));
205 hlb->GetYaxis()->SetRangeUser(0, Double_t(lbn_range));
206 }
207 TH1I* hlbstat(0);
208 for (std::vector<std::string>::iterator it = lbnDirs.begin(); it != lbnDirs.end(); ++it, ++k) {
209 if (k > nbinx) continue;
210 std::string histLBstat = run_dir + std::string("/") + *it + std::string("/HLT") + (*itFex) + "/lbnstatus/EF_MET_status";
211 if (!f->Get(histLBstat.c_str())) {
212 //if (dbgLevel) std::cerr << "--> HLTMETDQFlagSummary: no histogram " << histLBstat << std::endl;
213 continue;
214 }
215 hlbstat = (TH1I*) (f->Get(histLBstat.c_str()));
216 int flag = HLTMETGetStatusPerBin(hlbstat, 17, 29, 32, 32);
217 //std::cout << "--> HLTMETDQFlagSummary: " << *it << "\t: flag = " << flag << std::endl;
218 TString label = *it;
219 label.ReplaceAll("lowStat_", "");
220 hlb->GetXaxis()->SetBinLabel(k, label.Data());
221 hlb->GetYaxis()->SetBinLabel(k, label.Data());
222 hlb->SetBinContent(k, k, flag);
223 } // end for
224 hlb->Write("", TObject::kOverwrite);
225 } else {
226 std::cerr << "--> HLTMETDQFlagSummary: lowStat_* directories: found none " << std::endl;
227 }
228 } // end loop over fexs
229 } // end while loop over all keys
230 } // end method MonitoringFile::HLTMETDQFlagSummary
231
232 //---------------------------------------------------------------------------------------------------
233 // Method to get lowStat_* directories
234 //---------------------------------------------------------------------------------------------------
235 size_t MonitoringFile::HLTMETGetDQLBNRange(TDirectory*& run_dir, std::vector<std::string>& lbnDirs) {
236 if (!run_dir) return 0;
237
238 lbnDirs.clear();
239
240 //bool dbgLevel = false;
241 //if (dbgLevel) std::cout << "--> HLTMETGetDQLBNRange: lowStat_* directories: ";
242 run_dir->cd();
243
244 TIter next_lbn(run_dir->GetListOfKeys());
245 TKey* key_lbn(0);
246 std::string lbn_dir = "";
247 // loop over all objects
248 while ((key_lbn = dynamic_cast<TKey*> (next_lbn())) != 0) {
249 if (!key_lbn->IsFolder()) continue;
250 lbn_dir = key_lbn->GetName();
251 if (lbn_dir.find("lowStat_") == std::string::npos) {
252 continue;
253 }
254 lbnDirs.push_back(lbn_dir);
255 //std::cout << "dir = " << lbn_dir << "\tuuid = " << uuid.AsString() << "\ttime = " << dtme.GetTime() <<
256 // std::endl;
257 }
258 unsigned int nLBNDirs = lbnDirs.size();
259 //if (dbgLevel) std::cout << "found " << nLBNDirs << std::endl;
260 return nLBNDirs;
261 } // end method MonitoringFile::HLTMETGetDQLBNRange
262
263 //---------------------------------------------------------------------------------------------------
264 // Method to Get Status flags per lbn block
265 //---------------------------------------------------------------------------------------------------
266 int MonitoringFile::HLTMETGetStatusPerBin(TH1I*& hist, int yellmin, int yellmax, int redmin, int redmax) {
267 if (!hist) return 0;
268
269 std::string hname = hist->GetName();
270 int flag = 10; // 1 == GREEN, 2 == YELLOW, 3 == RED
271 float epsilon = 1.e-3;
272 int nbins = (int) hist->GetNbinsX();
273 yellmin = (yellmin > 0 && yellmin <= nbins) ? yellmin : 1;
274 yellmax = (yellmax > 0 && yellmax <= nbins) ? yellmax : nbins;
275 redmin = (redmin > 0 && redmin <= nbins) ? redmin : 1;
276 redmax = (redmax > 0 && redmax <= nbins) ? redmax : nbins;
277 for (int k = 1; k <= nbins; k++) {
278 float theVal = hist->GetBinContent(k);
279 // at least one bin > 0
280 if (theVal > epsilon) {
281 if (k >= yellmin && k <= yellmax) flag = 100;
282 if (k >= redmin && k <= redmax) flag = 1000;
283 } // end if
284 } // end for
285 return flag;
286 } // end method MonitoringFile::HLTMETGetStatusPerBin
287} // end namespace
static int HLTMETGetStatusPerBin(TH1I *&hist, int ymin, int ymax, int rmin, int rmax)
static size_t HLTMETGetDQLBNRange(TDirectory *&run_dir, std::vector< std::string > &lbnDirs)
static void HLTMETDQFlagSummary(TFile *f, std::string &run_dir)
std::string label(const std::string &format, int i)
Definition label.h:19