ATLAS Offline Software
MonitoringFile_HLTTauPostProcess.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3  */
4 
5 // **********************************************************************
6 // HLT Tau post processing jobs
7 // 01/03/2011: update by ccuenca
8 // **********************************************************************
9 
11 
12 #include <cmath>
13 #include <vector>
14 
15 #include <TCanvas.h>
16 #include <TF1.h>
17 #include <TFile.h>
18 #include <TH1.h>
19 #include <TH2.h>
20 #include <TKey.h>
21 #include <TMath.h>
22 #include <TProfile.h>
23 
24 namespace dqutils {
25  //--------------------------------------------------------------------------------
26  // main process
27  //--------------------------------------------------------------------------------
28  void MonitoringFile::HLTTauPostProcess(const std::string& inFilename, bool /* isIncremental */) {
29  //std::cout << "--> HLTTauPostProcess: Begin HLT Tau post-processing" << std::endl;
30 
31  //open root file
32  TFile* f = TFile::Open(inFilename.c_str(), "UPDATE");
33 
34  //check files are loaded.
35  if (f == 0 || !f->IsOpen()) {
36  //std::cerr << "--> HLTTauPostProcess: Input file not opened" << std::endl;
37  delete f;
38  return;
39  }
40 
41  //zombie files?
42  if (f->IsZombie()) {
43  //std::cerr << "--> HLTTauPostProcess: Input file " << inFilename << " cannot be opened. " << std::endl;
44  delete f;
45  return;
46  }
47 
48  //check file size is not too small.
49  if (f->GetSize() < 1000.) {
50  //std::cerr << "--> HLTTauPostProcess: Input file empty" << std::endl;
51  f->Close();
52  delete f;
53  return;
54  }
55 
56 
57  //making sure we are in the right directory
58  // f->cd("/");
59  // f->pwd();
60  // f->ls();
61 
62  //build iterator
63  TIter next_run(f->GetListOfKeys());
64  TKey* key_run(0);
65 
66  //loop over keys in root directory
67  while ((key_run = dynamic_cast<TKey*>(next_run())) != 0) {
68  TObject* obj_run = key_run->ReadObj();
69  TDirectory* tdir_run = dynamic_cast<TDirectory*>(obj_run);
70 
71  //check we have a valid pointer
72  if (tdir_run == 0) {
73  delete obj_run;
74  continue;
75  }
76 
77  //check name
78  std::string runDirName(tdir_run->GetName());
79  //std::cout<<"Run_directory: "<<runDirName<<std::endl;
80 
81  if (runDirName.find("run") == std::string::npos) {
82  delete obj_run;
83  continue;
84  }
85 
86  //find tauMon dir
87  TString taumonDirName = runDirName + "/HLT/TauMon";
88  TDirectory* taumonDir(0);
89  if (!(taumonDir = f->GetDirectory(taumonDirName))) {
90  //std::cerr << "--> HLTTauPostProcess: directory " << taumonDirName << " not found." << std::endl;
91  return;
92  }
93 
94 
95  std::vector<TString> varName;
96  varName.push_back("Et");
97  varName.push_back("Eta");
98  varName.push_back("Phi");
99 
100  std::vector<TString> lvlName;
101  lvlName.push_back("L1");
102  lvlName.push_back("L2");
103  lvlName.push_back("EF");
104 
105  std::vector< std::pair< int, int > > ratioIndex;
106  ratioIndex.push_back(std::make_pair(2, 0)); // EF vs L1
107  ratioIndex.push_back(std::make_pair(1, 0)); // L2 vs L1
108  ratioIndex.push_back(std::make_pair(2, 1)); // EF vs L2
109 
110 
111 
112  std::vector<TString> varName0;
113  varName0.push_back("Pt");
114  varName0.push_back("Eta");
115  varName0.push_back("Phi");
116  varName0.push_back("Nvtx");
117 
118  std::vector<TString> lvlNameO;
119  lvlNameO.push_back("L1");
120  lvlNameO.push_back("L2");
121  lvlNameO.push_back("EF");
122  lvlNameO.push_back("");//offline
123 
124  std::vector< std::pair< int, int > > ratioIndexO;
125  ratioIndexO.push_back(std::make_pair(0, 3)); // L1 vs off
126  ratioIndexO.push_back(std::make_pair(1, 3)); // L2 vs off
127  ratioIndexO.push_back(std::make_pair(2, 3)); // EF vs off
128 
129 
130  //loop over directories
131  TIter next_trig(taumonDir->GetListOfKeys());
132  TKey* key_trig(0);
133  while ((key_trig = dynamic_cast<TKey*>(next_trig())) != 0) {
134  TObject* obj_trig = key_trig->ReadObj();
135  TDirectory* dir_trig = dynamic_cast<TDirectory*>(obj_trig);
136  if (!dir_trig) continue;
137  //std::cout<<"--> HLTTauPostProcess: calling functions for " << dir_trig->GetName() << endl;
138  // "<<runDirName<<" and trigger item="<<dir_trig->GetName()<<endl;
139 
140  HLTTauPostProcess(f, dir_trig, "/RelativeEfficiency", "/RelativeEfficiency/Efficiency",
141  lvlName, varName, ratioIndex, 1);//name style : 1 (relative eff)
142  HLTTauPostProcess(f, dir_trig, "/OfflineRatio", "/OfflineRatio/Ratio",
143  lvlNameO, varName0, ratioIndexO, 2);//name style : 2 (offline eff)
144  HLTTauPostProcess(f, dir_trig, "/OfflineRatio/BDTMedium", "/OfflineRatio/Ratio/BDTMedium",
145  lvlNameO, varName0, ratioIndexO, 2, "BDTMedium");//name style : 2 (offline eff)
146  }
147  }
148 
149  f->Close();
150  delete f;
151  //std::cout << "--> HLTTauPostProcess: finished HLT Tau post-processing"<<std::endl;
152  }
153 
154  void MonitoringFile::HLTTauPostProcess(TFile* f, TDirectory* dir,
155  TString pathApp, TString pathAppEff,
156  const std::vector<TString>& lvlN, const std::vector<TString>& varN,
157  const std::vector< std::pair< int, int > >& ratioIndex, int nameStyle,
158  TString nameApp) {
159  std::string path = getPath(dir);
160 
161  //cout<<"HLTTauPostProcess in "<< path
162 // <<" pathApp="<<pathApp
163 // <<" pathAppEff="<<pathAppEff
164 // <<" nameStyle="<<nameStyle<<endl;
165 
166  TString basePath = path + pathApp + "/";
167  if (f->cd(basePath.Data()) == 0) {
168  //cout<<"basePath isn't there!"<<endl;
169  return;
170  }
171 
172  //TH1F* hRoI[lvlN.size()][varN.size()];
173  TH1F* hRoI[100][100];
174  for (unsigned int iLvl = 0; iLvl < lvlN.size(); iLvl++) {
175  for (unsigned int iVar = 0; iVar < varN.size(); iVar++) {
176  TString hName;
177  if (nameStyle == 1) hName = basePath + "h" + lvlN[iLvl] + "RoI" + varN[iVar] + (iLvl == 0 ? "Denom" : "Num") + nameApp;
178  if (nameStyle == 2) hName = basePath + "hTau" + varN[iVar] + lvlN[iLvl] + nameApp;
179  if (!CheckHistogram(f, hName.Data())) {
180  //cout<<" histo "<<hName<<" is not in f "<<f->GetName()<<endl;
181  return;
182  }
183  hRoI[iLvl][iVar] = (TH1F*) (f->Get(hName.Data()))->Clone();
184  }
185  }
186 
187 
188  basePath += path + pathAppEff + "/";
189  f->cd(basePath.Data());
190  TH1F* hEff[100][100];
191  for (unsigned int iVar = 0; iVar < varN.size(); iVar++) {
192  for (unsigned int iRatio = 0; iRatio < ratioIndex.size(); iRatio++) {
193  TString hName;
194  if (nameStyle == 1)
195  hName = basePath + "h" +
196  lvlN[ ratioIndex[ iRatio ].first ] + "vs" + lvlN[ ratioIndex[ iRatio ].second ] +
197  varN[iVar] + "Efficiency" + nameApp;
198 
199  if (nameStyle == 2)
200  hName = basePath + "h" +
201  lvlN[ ratioIndex[ iRatio ].first ] + "OfflineRatio" +
202  varN[iVar] + nameApp;
203 
204  if (!CheckHistogram(f, hName.Data())) return;
205 
206  hEff[iRatio][iVar] = (TH1F*) (f->Get(hName.Data()));
207  hEff[iRatio][iVar]->Divide(hRoI[ ratioIndex[ iRatio ].first ][iVar],
208  hRoI[ ratioIndex[ iRatio ].second ][iVar],
209  1.0, 1.0, "b");
210 
211  hEff[iRatio][iVar]->Write("", TObject::kOverwrite);
212  }
213  }
214 
215  f->Write();
216  }
217 }
python.SystemOfUnits.second
int second
Definition: SystemOfUnits.py:120
python.CaloRecoConfig.f
f
Definition: CaloRecoConfig.py:127
athena.path
path
python interpreter configuration --------------------------------------—
Definition: athena.py:126
dqutils::MonitoringFile::CheckHistogram
static bool CheckHistogram(TFile *f, const char *HistoName)
PixelAthClusterMonAlgCfg.varName
string varName
end cluster ToT and charge
Definition: PixelAthClusterMonAlgCfg.py:117
dqutils::MonitoringFile::getPath
static std::string getPath(TDirectory *dir)
dqutils
Definition: CoolMdt.h:76
dqutils::MonitoringFile::HLTTauPostProcess
static void HLTTauPostProcess(const std::string &inFilename, bool isIncremental=false)
Definition: MonitoringFile_HLTTauPostProcess.cxx:28
beamspotman.dir
string dir
Definition: beamspotman.py:623
MonitoringFile.h
TH1F
Definition: rootspy.cxx:320
DeMoScan.first
bool first
Definition: DeMoScan.py:534