ATLAS Offline Software
Loading...
Searching...
No Matches
MonitoringFile_HLTCaloAveragePtPhiEtaMaps.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/* HLTCalo Post Processing Method: to produce eta phi map of the average hit energy in each bin, made
6 * from two 2-d histograms, one with total energy in each bin the
7 * other with number of hits per bin
8 * Author : Gareth Brown <gabrown@CERN.CH> based on HLTMET Post Processing
9 * Date : Oct 2011
10 */
11
13
14#include <iostream>
15#include <cstdlib>
16#include <vector>
17#include <string>
18
19#include "TH2F.h"
20#include "TFile.h"
21#include "TKey.h"
22
23namespace dqutils {
24 void MonitoringFile::HLTCaloAveragePtPhiEtaMaps(TFile* f, TString& run_dir) {
25 int debugLevel = MonitoringFile::getDebugLevel();
26
27 if (debugLevel > 1) std::cout << "--> HLTCaloAveragePtPhiEtaMaps: <Quantity(eta,phi)> = Quantity(eta,phi)/N(eta,phi) " << std::endl;
28
29 f->cd("/");
30 TIter next_run(f->GetListOfKeys());
31 TKey* key_run(0);
32 while ((key_run = dynamic_cast<TKey*> (next_run())) != 0) {
33 if (!key_run->IsFolder()) continue;
34 run_dir = key_run->GetName();
35 if (!run_dir.Contains("run")) {
36 continue;
37 }
38
39 std::string run_dir2 = run_dir.Data();
40
41 // all merged root files have the structure "rootfile.root:/run_NNNNNNN"
42 // use that to extract run number
43 //int run_number = atoi( (run_dir2.substr(4, run_dir2.size()-4 )).c_str() );
44 //run_number=run_number;
45
46 // begin HLTMET
47 // note 1: prefix all dirs and hists with '/'
48 // note 2: missing dir => return
49 // note 3: missing hist => continue
50 TString hlt_top = run_dir + "/HLT"; // toplevel
51
52 std::vector<TString> calo_fold;
53 std::vector<std::pair<TString, TString> > hist_numr;
54
55 // The folders that we want to get the hists from
56 calo_fold.push_back("/CaloMon");
57 calo_fold.push_back("/CaloMonL2");
58
59 // check if the folders are in hlt
60 for (std::vector<TString>::iterator it = calo_fold.begin(); it != calo_fold.end(); ++it) {
61 TString theHistDir = hlt_top + *it;
62 TDirectory* dir = f->GetDirectory(theHistDir);
63 if (!dir) {
64 if (debugLevel > 0) std::cerr << "--> HLTCaloAveragePtPhiEtaMaps: directory " << theHistDir << " not found" << std::endl;
65 return;
66 }
67 }
68
69 // pairs of Num and Dem
70 hist_numr.push_back(std::make_pair<TString, TString>("/EnergyAccetaphiLAr", "/HitAccetaphiLAr"));
71 hist_numr.push_back(std::make_pair<TString, TString>("/EnergyAccetaphiTile", "/HitAccetaphiTile"));
72
73
74 // we have all dirs, get the component histograms
75 for (std::vector<TString>::iterator itFex = calo_fold.begin(); itFex != calo_fold.end(); ++itFex) {
76 for (std::vector<std::pair<TString, TString> >::iterator itNum = hist_numr.begin(); itNum != hist_numr.end(); ++itNum) {
77 TH2F* hnum(0), *hden(0);
78
79 // prepend histogram name with path and append with suffix [_00 .., _24 for each component]
80 TString thePath = hlt_top + (*itFex);
81 TString numPath = thePath + (*itNum).first;
82 TString denPath = thePath + (*itNum).second;
83
84 // get histograms
85 hnum = (TH2F*) (f->Get(numPath));
86 hden = (TH2F*) (f->Get(denPath));
87 // test if histograms are present
88 if (!hnum) {
89 if (debugLevel > 0) std::cerr << "--> HLTCaloAveragePtPhiEtaMaps: no histogram " << numPath << std::endl;
90 continue;
91 }
92 if (!hden) {
93 if (debugLevel > 0) std::cerr << "--> HLTCaloAveragePtPhiEtaMaps: no histogram " << denPath << std::endl;
94 continue;
95 }
96
97
98 // get directory of histograms
99 TDirectory* dir = f->GetDirectory(thePath);
100
101 // these are disabled, because we have to worry about updating metadata
102 // clone numerator histogram in the same directory; prepend with "avg_"
103 // use numerator to do the job
104
105 // divide num by den to get average quantity.
106 hnum->Divide(hnum, hden);
107
108 // fix histogram titles
109 TString title = hnum->GetTitle();
110 hnum->SetTitle("Average Transverse Energy per eta/phi bin");
111
112 dir->cd();
113 hnum->Write("", TObject::kOverwrite);
114 } // done looping over quantities
115 } // done looping over directories
116 } // end while loop over all keys
117 } // end method MonitoringFile::HLTCaloAveragePtPhiEtaMaps
118} // end namespace
static void HLTCaloAveragePtPhiEtaMaps(TFile *f, TString &run_dir)
static int getDebugLevel()