ATLAS Offline Software
postProcessIDPVMHistos.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2021 CERN for the benefit of the ATLAS collaboration
3 */
4 
17 
18 #include <iostream>
19 
20 #include "TFile.h"
21 #include "TSystem.h"
22 #include "TH1.h"
23 #include "TH2.h"
24 #include "TObject.h"
25 
26 using namespace std;
27 
28 
29 bool file_exists(const string & p_name) {
30  return !gSystem->AccessPathName(p_name.c_str(), kFileExists);
31 }
32 
33 // check if the name of an object matches what we expect from a resolution helper
34 bool isResolutionHelper(TObject* entry){
35  const std::string objName{entry->GetName()};
36  return ((objName.find("resHelper") == 0 || objName.find("pullHelper") == 0) && dynamic_cast<TH1*>(entry));
37 }
38 
39 // get the x axis observable and resolution (y axis) from the name of a 2D histo.
40 // Relies on the conventions within IDPVM
41 std::pair<std::string, std::string> getObservableAndReso(const TObject* resHelper){
42  const std::string name{resHelper->GetName()};
43  const std::string keyWord {"Helper_"};
44  const size_t offset = keyWord.size();
45  auto start = name.find(keyWord)+offset;
46  auto sep = name.find('_',start);
47  return {name.substr(start, sep-start), name.substr(sep+1)};
48 
49 }
50 
51 // get the resolution type (pull or res) - taken from the prefix of the 2D histo name
52 std::string getResoType(const TObject* resHelper){
53  std::string aux{resHelper->GetName()};
54  return aux.substr(0, aux.find("Helper"));
55 }
56 
57 // clone an existing histogram of a known name
58 TH1* cloneExisting(const std::string & name){
59  auto *h = gDirectory->Get(name.c_str());
60  if (!h){
61  std::cerr << "Could not find existing histogram "<<name<<" - will not postprocess "<<std::endl;
62  return nullptr;
63  }
64  auto *ret = dynamic_cast<TH1*>(h->Clone(name.c_str()));
65  if (!ret){
66  std::cerr << "Found an existing object "<<name<<", but it is not a histogram ("<<h->IsA()->GetName()<<") - will not postprocess "<<std::endl;
67  }
68  return ret; // will also catch ret == nullptr
69 }
70 
71 // get the names of the 1D histograms following IDPVM conventions.
72 std::pair<std::string, std::string> getPullAndResoNames(const std::string & type){
73  if (type == "res"){
74  return {"resolution","resmean"};
75  }
76  else if (type == "pull"){
77  return {"pullwidth","pullmean"};
78  }
79  else {
80  std::cerr << " Not able to identify the histogram names for a resolution type "<<type<<" - supported are 'res' and 'pull'. "<<std::endl;
81  }
82  return {"",""};
83 }
84 
85 int postProcessHistos(TObject* resHelper, IDPVM::ResolutionHelper & theHelper){
86  // here we have to rely on the naming conventions of IDPVM to identify what we are looking at
87  auto vars = getObservableAndReso(resHelper);
88  auto type = getResoType(resHelper);
89  // cast to TH2
90  TH2* resHelper2D = dynamic_cast<TH2*>(resHelper);
91  if (!resHelper2D){
92  std::cerr <<"Unable to reduce the histogram "<<resHelper->GetName()<<" to a TH2 - this histo can not yet be postprocessed! " <<std::endl;
93  return 1;
94  }
95  const auto & oneDimNames = getPullAndResoNames(type);
96  // get the corresponding 1D histos by cloning the existing ones in the same folder
97  TH1* h_width = cloneExisting(oneDimNames.first+"_vs_"+vars.first+"_"+vars.second);
98  TH1* h_mean = cloneExisting(oneDimNames.second+"_vs_"+vars.first+"_"+vars.second);
99  // then call the resolution helper as done in "online" IDPVM
100  theHelper.makeResolutions(resHelper2D, h_width, h_mean);
101  // update our 1D histos
102  h_width->Write();
103  h_mean->Write();
104  // and we are done
105  return 0;
106 }
107 
108 // recursively parse a directory tree, post-processing any histos seen along the way
109 int postProcessDir(TDirectory* dir, IDPVM::ResolutionHelper & theHelper){
110 
111  int outcome = 0;
112  auto theCWD = gDirectory;
113  // walk through all keys in this directory
114  dir->cd();
115  auto *keys = dir->GetListOfKeys();
116  for (auto *const key : *keys){
117  TObject* gotIt = dir->Get(key->GetName());
118 
119  // if we encounter a directory, descend into it and repeat the process
120  TDirectory* theDir = dynamic_cast<TDirectory*>(gotIt);
121  if (theDir){
122  outcome |= postProcessDir(theDir, theHelper);
123  }
124 
125  // if we encounter a histogram that could be a resolution input, post-process it
126  if (isResolutionHelper(gotIt)){
127  outcome |= postProcessHistos(gotIt, theHelper);
128  }
129  }
130  theCWD->cd();
131  return outcome;
132 }
133 
134 // function driving the postprocessing for this file
135 int pproc_file(const std::string & p_infile) {
136 
137  IDPVM::ResolutionHelper theHelper;
138 
139  TFile* infile = TFile::Open(p_infile.c_str(),"UPDATE");
140  if (!infile ) {
141  std::cerr << "could not open input file "<<p_infile<<" for updating "<< std::endl;
142  return 1;
143  }
144 
145  int res = postProcessDir(gDirectory,theHelper); // recursively post-process the directory tree, starting from the root.
146  infile->Close();
147  return res;
148 }
149 
150 
151 int main(int argc, char* argv[]) {
152 
153  std::string infile{""};
155  if (argc == 2){
156  infile = argv[1];
157  }
162  else if (argc == 4 && std::string(argv[1]) == std::string{"-f"}){
164  gSystem->CopyFile(argv[3],argv[2],true);
166  infile = argv[2];
167  }
168  else {
169  std::cerr<<" Usage: postProcessIDPVMHistos <File to post-process>"<<std::endl;
170  std::cerr<< " where the file is typically obtained by hadding" << std::endl;
171  std::cerr<< " outputs of several independent IDPVM runs." << std::endl;
172  std::cerr<<" Alternative usage: postProcessIDPVMHistos -f <desired output file name> <File to post-process>"<<std::endl;
173  std::cerr<< " imitates a hadd-like signature for PhysVal merging." << std::endl;
174  return 1;
175  }
177  if (!file_exists(infile)) {
178  std::cerr << "Error: invalid input file: " << infile << std::endl;
179  return 1;
180  }
182  std::cout << " Post-processing file " << infile << "\n" << std::endl;
183  return pproc_file(infile);
184 }
isResolutionHelper
bool isResolutionHelper(TObject *entry)
Definition: postProcessIDPVMHistos.cxx:34
run.infile
string infile
Definition: run.py:13
IDPVM::ResolutionHelper::makeResolutions
void makeResolutions(const TH2 *h_input2D, TH1 *hwidth, TH1 *hmean, TH1 *hproj[], bool saveProjections, IDPVM::ResolutionHelper::methods theMethod=IDPVM::ResolutionHelper::iterRMS_convergence)
extract 1D resolution plots from a 2D "residual vs observable" histogram.
Definition: InnerDetector/InDetValidation/InDetPhysValMonitoring/src/ResolutionHelper.cxx:288
mergePhysValFiles.start
start
Definition: DataQuality/DataQualityUtils/scripts/mergePhysValFiles.py:14
postProcessHistos
int postProcessHistos(TObject *resHelper, IDPVM::ResolutionHelper &theHelper)
Definition: postProcessIDPVMHistos.cxx:85
cloneExisting
TH1 * cloneExisting(const std::string &name)
Definition: postProcessIDPVMHistos.cxx:58
pproc_file
int pproc_file(const std::string &p_infile)
Definition: postProcessIDPVMHistos.cxx:135
IDPVM::ResolutionHelper
Definition: InnerDetector/InDetValidation/InDetPhysValMonitoring/InDetPhysValMonitoring/ResolutionHelper.h:28
ResolutionHelper.h
getPullAndResoNames
std::pair< std::string, std::string > getPullAndResoNames(const std::string &type)
Definition: postProcessIDPVMHistos.cxx:72
ParseInputs.gDirectory
gDirectory
Definition: Final2012/ParseInputs.py:133
h
LArCellNtuple.argv
argv
Definition: LArCellNtuple.py:152
res
std::pair< std::vector< unsigned int >, bool > res
Definition: JetGroupProductTest.cxx:14
file_exists
bool file_exists(const string &p_name)
Definition: postProcessIDPVMHistos.cxx:29
DQHistogramMergeRegExp.argc
argc
Definition: DQHistogramMergeRegExp.py:20
main
int main(int argc, char *argv[])
Definition: postProcessIDPVMHistos.cxx:151
beamspotman.dir
string dir
Definition: beamspotman.py:623
GetAllXsec.entry
list entry
Definition: GetAllXsec.py:132
grepfile.sep
sep
Definition: grepfile.py:38
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:221
getObservableAndReso
std::pair< std::string, std::string > getObservableAndReso(const TObject *resHelper)
Definition: postProcessIDPVMHistos.cxx:41
python.CaloScaleNoiseConfig.type
type
Definition: CaloScaleNoiseConfig.py:78
convertTimingResiduals.offset
offset
Definition: convertTimingResiduals.py:71
python.Bindings.keys
keys
Definition: Control/AthenaPython/python/Bindings.py:798
getResoType
std::string getResoType(const TObject *resHelper)
Definition: postProcessIDPVMHistos.cxx:52
postProcessDir
int postProcessDir(TDirectory *dir, IDPVM::ResolutionHelper &theHelper)
Definition: postProcessIDPVMHistos.cxx:109
mapkey::key
key
Definition: TElectronEfficiencyCorrectionTool.cxx:37