ATLAS Offline Software
MonitoringFile_HLTMinBiasMonPostProcess.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3  */
4 
5 // **********************************************************************
6 // HLT MinBias post processing jobs
7 // 09/06/2015: update by Andrzej Zemla <azemla@cern.ch>
8 // **********************************************************************
9 
11 
12 #include <cmath>
13 #include <vector>
14 
15 #include <TF1.h>
16 #include <TFile.h>
17 #include <TH1.h>
18 #include <TKey.h>
19 #include <TMath.h>
20 #include <TROOT.h>
21 #include <TCanvas.h>
22 #include <TString.h>
23 
24 namespace dqutils {
25  //--------------------------------------------------------------------------------
26  // main process
27  //--------------------------------------------------------------------------------
28  void MonitoringFile::HLTMinBiasMonPostProcess(const std::string& inFilename, bool /* isIncremental */) {
29  // std::cout << "--> HLTMinBiasMonPostProcess: Begin HLT MinBias post-processing" << std::endl;
30 
31  //open root file
32  TFile* f = TFile::Open(inFilename.c_str(), "READ");
33 
34  //check files are loaded.
35  if (f == 0 || !f->IsOpen()) {
36  // std::cerr << "--> HLTMinBiasMonPostProcess: Input file not opened" << std::endl;
37  delete f;
38  return;
39  }
40 
41  //zombie files?
42  if (f->IsZombie()) {
43  // std::cerr << "--> HLTMinBiasMonPostProcess: Input file " << inFilename << " cannot be opened. " <<
44  // std::endl;
45  delete f;
46  return;
47  }
48 
49  //check file size is not too small.
50  if (f->GetSize() < 1000.) {
51  // std::cerr << "--> HLTMinBiasMonPostProcess: Input file empty" << std::endl;
52  f->Close();
53  delete f;
54  return;
55  }
56  //build iterator
57  TIter next_run(f->GetListOfKeys());
58  TKey* key_run(0);
59  TString minbiasmonDirName;
60 
61  std::vector< std::pair<TString, TString> >* v_targetNames = new std::vector< std::pair<TString, TString> >(0);
62 
63 
64  //loop over keys in root directory
65  while ((key_run = dynamic_cast<TKey*>(next_run())) != 0) {
66  TObject* obj_run = key_run->ReadObj();
67  TDirectory* tdir_run = dynamic_cast<TDirectory*>(obj_run);
68 
69  //check we have a valid pointer
70  if (tdir_run == 0) {
71  delete obj_run;
72  continue;
73  }
74 
75  //check name
76  std::string runDirName(tdir_run->GetName());
77  // std::cout<<"Run_directory: "<<runDirName<<std::endl;
78 
79  if (runDirName.find("run") == std::string::npos) {
80  delete obj_run;
81  continue;
82  }
83 
84  //find MinBiasMon dir
85  minbiasmonDirName = runDirName + "/HLT/MinBiasMon";
86  TDirectory* minbiasmonDir(0);
87  if (!(minbiasmonDir = f->GetDirectory(minbiasmonDirName))) {
88 // std::cerr << "--> HLTMinBiasMonPostProcess: directory " << minbiasmonDirName << " not found." << std::endl;
89  delete v_targetNames;
90  return;
91  }
93  }
94 
95  f->Close();
96  f = 0;
97 
98  f = TFile::Open(inFilename.c_str(), "UPDATE");
99 
100  const TString repPath = TString(f->GetPath());
101 
102  //first update global efficiency
103  f->cd(minbiasmonDirName);
104 
105  TH1F* h_triggEffic = dynamic_cast<TH1F*>(gDirectory->Get("TriggerEfficiencies"));
106  TH1F* h_triggEfficPass = dynamic_cast<TH1F*>(gDirectory->Get("TriggerEfficienciesPassed"));
107  TH1F* h_triggEfficAll = dynamic_cast<TH1F*>(gDirectory->Get("TriggerEfficienciesAll"));
108  if (!h_triggEffic or !h_triggEfficPass or !h_triggEfficAll) {
109  std::cerr << "Dynamic cast failed in MonitoringFile::HLTMinBiasMonPostProcess\n";
110  f->Close();
111  delete f;
112  delete v_targetNames;
113  return;
114  }
115 
116  h_triggEffic->Divide(h_triggEfficPass, h_triggEfficAll, 1., 1., "B");
117  h_triggEffic->GetYaxis()->SetRangeUser(0., 1.2);
118  h_triggEffic->Write("", TObject::kOverwrite);
119 
120  for (uint k = 0; k < v_targetNames->size(); k++) {
121  v_targetNames->at(k).first.ReplaceAll(repPath, "");
122 
123  f->cd(v_targetNames->at(k).first);
124 
125  TH1F* h_target = dynamic_cast<TH1F*>(gDirectory->Get(v_targetNames->at(k).second));
126  TH1F* h_num = dynamic_cast<TH1F*>(gDirectory->Get(v_targetNames->at(k).second + "Passed"));
127  TH1F* h_den = dynamic_cast<TH1F*>(gDirectory->Get(v_targetNames->at(k).second + "All"));
128 
129  if (h_target != 0 && h_num != 0 && h_den != 0) {
130  h_target->Divide(h_num, h_den, 1., 1., "B");
131  h_target->GetYaxis()->SetRangeUser(0., 1.2);
132  h_target->Write("", TObject::kOverwrite);
133  }
134  }
135  v_targetNames->clear();
136  delete v_targetNames;
137  v_targetNames = 0;
138 
139  f->Close();
140  delete f;
141  }
142 
143  void MonitoringFile::HLTMinBiasMonGetTargetHistos(TDirectory* source, std::vector< std::pair< TString,
144  TString > >& targetNames)
145  {
146  TKey* key;
147  TIter nextkey(source->GetListOfKeys());
148 
149  while ((key = (TKey*) nextkey())) {
150  const char* classname = key->GetClassName();
151  TClass* cl = gROOT->GetClass(classname);
152  if (!cl) continue;
153  if (cl->InheritsFrom(TDirectory::Class())) {
154  TString kname = key->GetName();
155 
156  source->cd(key->GetName());
157  TDirectory* nextdir = gDirectory;
158  HLTMinBiasMonGetTargetHistos(nextdir, targetNames);
159  }
160  if (cl->InheritsFrom(TH1F::Class())) {
161  TString cRatio = TString(key->GetName());
162  if (cRatio == "Efficiency" ||
163  cRatio == "Purity" ||
164  cRatio == "EfficienciesTrigger" ||
165  cRatio == "EfficiencyTracks" ||
166  cRatio ==
167  "TriggerPurities") targetNames.push_back(std::pair<TString, TString>(TString(gDirectory->GetPath()),
168  TString(key->GetName())));
169  }
170  }
171  }
172 }
python.CaloRecoConfig.f
f
Definition: CaloRecoConfig.py:127
uint
unsigned int uint
Definition: LArOFPhaseFill.cxx:20
dqutils::MonitoringFile::HLTMinBiasMonPostProcess
static void HLTMinBiasMonPostProcess(const std::string &inFileName, bool isIncremental=false)
Definition: MonitoringFile_HLTMinBiasMonPostProcess.cxx:28
ParseInputs.gDirectory
gDirectory
Definition: Final2012/ParseInputs.py:133
dqutils
Definition: CoolMdt.h:76
MonitoringFile.h
TH1F
Definition: rootspy.cxx:320
dqutils::MonitoringFile::HLTMinBiasMonGetTargetHistos
static void HLTMinBiasMonGetTargetHistos(TDirectory *source, std::vector< std::pair< TString, TString > > &targetNames)
Definition: MonitoringFile_HLTMinBiasMonPostProcess.cxx:143
dq_make_web_display.cl
cl
print [x.__class__ for x in toList(dqregion.getSubRegions()) ]
Definition: dq_make_web_display.py:26
fitman.k
k
Definition: fitman.py:528
mapkey::key
key
Definition: TElectronEfficiencyCorrectionTool.cxx:37