ATLAS Offline Software
MonitoringFile_HLTEgammaEfficiency.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3  */
4 
5 // **********************************************************************
6 // HLT Egamma post processing jobs
7 // 15/03/2012 initial file by Yan-Jie Schnellbach, based on HLT Tau
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  // function to calculate bootstrap efficiency w.r.t to offline
27  // ------------------------------------------------------------
28  void MonitoringFile::HLTEgammaEfficiencyOff(TFile* f, TDirectory* trig_dir,
29  TDirectory* off_dir, const TString& pathNum,
30  const TString& pathEff, const std::vector<TString>& varName) {
31  //get the relevant paths
32  std::string chain_path = getPath(trig_dir);
33  std::string offline_path = getPath(off_dir);
34 
35  //establish base path (to save typing "/PassedChain/" over and over)
36  TString basePath = chain_path + "/PassedChain/";
37  if (f->cd(basePath.Data()) == 0) {
38  std::cerr << "----> HLTEgammaPostProcess: basePath not found in "
39  << chain_path << ", skipping now!\n";
40  return;
41  }
42 
43  //create histogram matrix with pointers to our histograms in question
44  // 3x3 matrix, rows are variables (Et, Eta, Phi), colums the different
45  // histograms (numerator - matched offline, denominator - offline, eff)
46  TH1F* histo_matrix[3][3];
47 
48  //check variables to avoid segfaults!
49  if (varName.size() > 3) {
50  std::cerr << "----> HLTEgammaPostProcess: Too many variables "
51  << "used in " << chain_path + "/PassedChain/" + pathEff
52  << " for bootstrap efficiencies, skipping now!\n";
53  return;
54  }
55 
56 
57 
58  for (uint iVar = 0; iVar < varName.size(); ++iVar) {
59  //form names for all variables
60  TString numName = basePath + pathNum + "eg" + varName[iVar];
61  TString denName = offline_path + "/" + "eg" + varName[iVar];
62  TString effName = basePath + pathEff + "eg" + varName[iVar];
63 
64  //check histogram existence
65  if (!CheckHistogram(f, numName.Data())) {
66  std::cerr << "----> HLTEgammaPostProcess: Histogram "
67  << numName << " does not exist in file "
68  << f->GetName() << " - skipping!\n";
69  return;
70  }
71  if (!CheckHistogram(f, denName.Data())) {
72  std::cerr << "----> HLTEgammaPostProcess: Histogram "
73  << denName << " does not exist in file "
74  << f->GetName() << " - skipping!\n";
75  return;
76  }
77  if (!CheckHistogram(f, effName.Data())) {
78  std::cerr << "----> HLTEgammaPostProcess: Histogram "
79  << effName << " does not exist in file "
80  << f->GetName() << " - skipping!\n";
81  return;
82  }
83 
84  //debug outs
85  //std::cout << "---->HLTEgammaPostProcess: Taking " << numName << " and divide by "
86  // << denName << " to get " << effName << std::endl;
87 
88  //now get histograms
89  histo_matrix[iVar][0] = (TH1F*) (f->Get(numName.Data()));
90  histo_matrix[iVar][1] = (TH1F*) (f->Get(denName.Data()));
91  histo_matrix[iVar][2] = (TH1F*) (f->Get(effName.Data()));
92 
93  //hop into efficiency path
94  TString writePath = basePath + pathEff;
95  f->cd(writePath.Data());
96 
97  //create efficiency histogram by division
98  histo_matrix[iVar][2]->Divide(histo_matrix[iVar][0], histo_matrix[iVar][1],
99  1.0, 1.0, "B");
100  histo_matrix[iVar][2]->Scale(100.);
101  if (histo_matrix[iVar][2]->Write("",
102  TObject::kOverwrite) ==
103  0) std::cerr << "----> HLTEgammaPostProcess: Efficiency histogram "
104  << effName <<
105  " not written! Histogram will be buggy!\n";
106  }
107 
108  //save results now
109  f->Write();
110  }
111 
112  // --------------------------------------------------
113  // function called to calculate relative efficiencies
114  // ---------------------------------------------------
115 
116  void MonitoringFile::HLTEgammaEfficiencyRel(TFile* f, TDirectory* trig_dir,
117  const TString& pathPre, const TString& pathRej,
118  const std::vector<TString>& objStage,
119  const std::vector<TString>& varName) {
120  std::string chain_path = getPath(trig_dir);
121 
122  // std::cout << "----> HLTEgammaPostProcess in " << chain_path
123  // << ", pathPre=" << chain_path + "/PassedChain/" + pathPre
124  // << ", pathRej=" << chain_path + "/PassedChain/" + pathRej
125 
126  //establish base path (to save typing)
127  TString basePath = chain_path + "/PassedChain/";
128  if (f->cd(basePath.Data()) == 0) {
129  std::cerr << "----> HLTEgammaPostProcess: basePath not found in "
130  << chain_path << ", skipping now!\n";
131  return;
132  }
133 
134  //create histogram matrix with pointers to our histograms in question
135  // 5x3x3 tensor, 1st index is the stage 2nd index are variables
136  // (Et, Eta, Phi), 3rd index the different histograms
137  // (numerator - n-1 stage, denominator - n stage)
138  TH1F* histo_tensor[5][3][3];
139 
140  //now check size of our stages and variables to avoid segfaults!
141  if (objStage.size() > 5) {
142  std::cerr << "----> HLTEgammaPostProcess: Too many chain stages "
143  << "found in " << chain_path << " for relative efficiencies, bailing it!\n";
144  return;
145  }
146  if (varName.size() > 3) {
147  std::cerr << "----> HLTEgammaPostProcess: Too many variables "
148  << "used in " << chain_path << " for relative efficiencies, bailing out!\n";
149  return;
150  }
151 
152  //iterate through variables (start at 1, as there is no we can't
153  //calculate rejection efficiency w.r.t n-1 with n = 0 -> segfault!)
154  for (uint iVar = 0; iVar < varName.size(); ++iVar) {
155  //iterate through stages
156  for (uint iStg = 1; iStg < objStage.size(); ++iStg) {
157  //form names for all variables
158  TString preName = basePath + pathPre + objStage[iStg - 1] + varName[iVar];
159  TString rejName = basePath + pathPre + objStage[iStg] + varName[iVar];
160  TString effName = basePath + pathRej + objStage[iStg] + varName[iVar];
161  //check histogram existence
162  if (!CheckHistogram(f, rejName.Data())) {
163  std::cerr << "----> HLTEgammaPostProcess: Histogram "
164  << rejName << " does not exist in file "
165  << f->GetName() << " - skipping it!\n";
166  return;
167  }
168  if (!CheckHistogram(f, preName.Data())) {
169  std::cerr << "----> HLTEgammaPostProcess: Histogram "
170  << preName << " does not exist in file "
171  << f->GetName() << " - skipping it!\n";
172  return;
173  }
174  if (!CheckHistogram(f, effName.Data())) {
175  std::cerr << "----> HLTEgammaPostProcess: Histogram "
176  << effName << " does not exist in file "
177  << f->GetName() << " - skipping it!\n";
178  return;
179  }
180 
181  //debug outs
182  //std::cout << "---->HLTEgammaPostProcess: Taking " << rejName << " and divide by "
183  // << preName << " to get " << effName << std::endl;
184 
185  //now get our histograms
186  histo_tensor[iStg][iVar][0] = (TH1F*) (f->Get(rejName.Data()));
187  histo_tensor[iStg][iVar][1] = (TH1F*) (f->Get(preName.Data()));
188  histo_tensor[iStg][iVar][2] = (TH1F*) (f->Get(effName.Data()));
189 
190  //hop into rejection path
191  TString writePath = basePath + pathRej;
192  f->cd(writePath.Data());
193 
194  //create rejection efficiency histograms by division
195  histo_tensor[iStg][iVar][2]->Divide(histo_tensor[iStg][iVar][0], histo_tensor[iStg][iVar][1],
196  1.0, 1.0, "B");
197  histo_tensor[iStg][iVar][2]->Scale(100.);
198  if (histo_tensor[iStg][iVar][2]->Write("", TObject::kOverwrite) == 0) std::cerr
199  << "---->HLTEgammaPostProcess: Relative efficiency histogram "
200  << effName << " not written! Histograms will be buggy!\n";
201  }
202  }
203 
204  //save results now
205  f->Write();
206  }
207 }
dqutils::MonitoringFile::CheckHistogram
static bool CheckHistogram(TFile *f, const char *HistoName)
dqutils::MonitoringFile::HLTEgammaEfficiencyOff
static void HLTEgammaEfficiencyOff(TFile *f, TDirectory *trig_dir, TDirectory *off_dir, const TString &pathNum, const TString &pathEff, const std::vector< TString > &varName)
Definition: MonitoringFile_HLTEgammaEfficiency.cxx:28
uint
unsigned int uint
Definition: LArOFPhaseFill.cxx:20
PixelAthClusterMonAlgCfg.varName
string varName
end cluster ToT and charge
Definition: PixelAthClusterMonAlgCfg.py:125
dqutils::MonitoringFile::HLTEgammaEfficiencyRel
static void HLTEgammaEfficiencyRel(TFile *f, TDirectory *trig_dir, const TString &pathPre, const TString &pathRej, const std::vector< TString > &objStage, const std::vector< TString > &varName)
Definition: MonitoringFile_HLTEgammaEfficiency.cxx:116
dqutils::MonitoringFile::getPath
static std::string getPath(TDirectory *dir)
hist_file_dump.f
f
Definition: hist_file_dump.py:135
dqutils
Definition: CoolMdt.h:76
MonitoringFile.h
python.TrigEgammaMonitorHelper.TH1F
def TH1F(name, title, nxbins, bins_par2, bins_par3=None, path='', **kwargs)
Definition: TrigEgammaMonitorHelper.py:24