ATLAS Offline Software
Loading...
Searching...
No Matches
MonitoringFile_HLTMETAveragePhivsEtaMaps.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4
5/* HLTMET Post Processing Method: Peform bin-wise division for each Phi.vs.Eta maps of MET, SumEt, SumE
6 * with Phi.vs.Eta map of N to compute <MET>(eta,phi), <SumEt>(eta,phi)
7 * and <SumE>(eta,phi)
8 * Author : Venkatesh Kaushik <venkat.kaushik@cern.ch>
9 * Date : Feb 2010
10 */
11
13
14#include <iostream>
15#include <iomanip>
16#include <algorithm>
17#include <fstream>
18#include <cmath>
19#include <cstdlib>
20#include <sstream>
21#include <vector>
22#include <utility>
23#include <map>
24#include <string>
25
26#include "TH1F.h"
27#include "TH2F.h"
28#include "TFile.h"
29#include "TClass.h"
30#include "TKey.h"
31#include "TMath.h"
32#include "TF1.h"
33#include "TTree.h"
34#include "TBranch.h"
35
36namespace dqutils {
37 void
38 MonitoringFile::HLTMETAveragePhivsEtaMaps(TFile* f, TString& run_dir) {
39 //bool dbgLevel = false;
40
41 //if (dbgLevel) std::cout << "--> HLTMETAveragePhivsEtaMaps: <Quantity(eta,phi)> = Quantity(eta,phi)/N(eta,phi) " << std::endl;
42
43 f->cd("/");
44 TIter next_run(f->GetListOfKeys());
45 TKey* key_run(0);
46 while ((key_run = dynamic_cast<TKey*> (next_run())) != 0) {
47 if (!key_run->IsFolder()) continue;
48 run_dir = key_run->GetName();
49 if (!run_dir.Contains("run")) {
50 continue;
51 }
52
53
54 std::string run_dir2 = run_dir.Data();
55
56 // all merged root files have the structure "rootfile.root:/run_NNNNNNN"
57 // use that to extract run number
58 //int run_number = atoi( (run_dir2.substr(4, run_dir2.size()-4 )).c_str() );
59 //run_number=run_number;
60
61 // begin HLTMET
62 // note 1: prefix all dirs and hists with '/'
63 // note 2: missing dir => return
64 // note 3: missing hist => continue
65 TString hlt_top = run_dir + "/HLT"; // toplevel
66 TString met_efdir = "/EFMissingET_Fex"; // EF dir
67
68 std::vector<TString> met_fexs, hist_numr;
69
70 // expect the following fex dirs
71 met_fexs.push_back("/METMon");
72 //met_fexs.push_back("/METMon_FEB");
73 //met_fexs.push_back("/METMon_allCells");
74
75 // check if fex dirs are in hlt
76 for (std::vector<TString>::iterator it = met_fexs.begin(); it != met_fexs.end(); ++it) {
77 TString theHistDir = hlt_top + *it;
78 TDirectory* dir = f->GetDirectory(theHistDir);
79 if (!dir) {
80 std::cerr << "--> HLTMETAveragePhivsEtaMaps: directory " << theHistDir << " not found" << std::endl;
81 return;
82 }
83 // expect the EF dir inside each fex dir
84 theHistDir += met_efdir;
85 dir = f->GetDirectory(theHistDir);
86 if (!dir) {
87 std::cerr << "--> HLTMETAveragePhivsEtaMaps: directory " << theHistDir << " not found" << std::endl;
88 return;
89 }
90 }
91
92 // MEt, SumEt, SumE of components (numerator)
93 // compSumEt_lin_EtaPhi_00 etc..
94 hist_numr.push_back("/compEt_lin_");
95 hist_numr.push_back("/compSumEt_lin_");
96 hist_numr.push_back("/compSumE_lin_"); // need to plot SumE on linear scale (todo)
97
98 // components N (denominator)
99 TString hist_denr = "/compN_";
100
101 // type (eta,phi map)
102 TString hist_suffix = "EtaPhi_"; // phi vs. eta map
103
104 // each component a 2d hist. get all components
105 unsigned int comp_num = 25; // 25 components
106
107 // we have all dirs, get the component histograms
108 for (std::vector<TString>::iterator itFex = met_fexs.begin(); itFex != met_fexs.end(); ++itFex) {
109 for (std::vector<TString>::iterator itNum = hist_numr.begin(); itNum != hist_numr.end(); ++itNum) {
110 for (unsigned int icomp = 0; icomp < comp_num; icomp++) {
111 TH2F* hnum(0), *hden(0);
112
113 // prepend histogram name with path and append with suffix [_00 .., _24 for each component]
114 TString thePath = hlt_top + (*itFex) + met_efdir;
115 TString numHist = (*itNum) + hist_suffix + TString(Form("%02u", icomp));
116 TString denHist = hist_denr + hist_suffix + TString(Form("%02u", icomp));
117 TString numPath = thePath + numHist;
118 TString denPath = thePath + denHist;
119
120 // test if histograms are present
121 if (!f->Get(numPath)) {
122 //if (dbgLevel) std::cerr << "--> HLTMETAveragePhivsEtaMaps: no histogram " << numPath << std::endl;
123 continue;
124 }
125 if (!f->Get(denPath)) {
126 //if (dbgLevel) std::cerr << "--> HLTMETAveragePhivsEtaMaps: no histogram " << denPath << std::endl;
127 continue;
128 }
129
130 // get histograms
131 hnum = (TH2F*) (f->Get(numPath));
132 hden = (TH2F*) (f->Get(denPath));
133
134 // get directory of histograms
135 TDirectory* dir = f->GetDirectory(thePath);
136
137 // these are disabled, because we have to worry about updating metadata
138 // clone numerator histogram in the same directory; prepend with "avg_"
139 // use numerator to do the job
140 //TString avgHist = TString("avg_") + (*itNum) + hist_suffix + TString(Form("%02u",icomp));
141 //havg = (TH2F *) (hnum->Clone(avgHist));
142 //havg->SetDirectory(dir);
143
144 // divide num by den to get average quantity.
145 hnum->Divide(hnum, hden);
146
147 // fix histogram titles
148 TString title = hnum->GetTitle();
149 title.ReplaceAll(": ", ": #LT");
150 title.ReplaceAll("(#eta", "#GT(#eta");
151 hnum->SetTitle(title);
152
153 //title = hden->GetTitle();
155 //title.ReplaceAll(" #phi VS #eta",": N(#eta, #phi)");
156 //hden->SetTitle(title);
157
158 //std::cout << numHist << " div " << denHist << std::endl;
159 //std::cout << hnum->GetZaxis()->GetTitle() << "\t";
160 //std::cout << hden->GetZaxis()->GetTitle() << std::endl;
161 //std::cout << hnum->GetTitle() << "\t";
162 //std::cout << hden->GetTitle() << std::endl;
163 //std::cout << "-------------------------------------------------------------" << std::endl;
164 //getchar();
165 // cd() into that directory and save the avg histogram in file
166 dir->cd();
167 hnum->Write("", TObject::kOverwrite);
168 } // done looping over components
169 } // done looping over quantities
170 } // done looping over directories
171 } // end while loop over all keys
172 } // end method MonitoringFile::HLTMETAveragePhivsEtaMaps
173} // end namespace
static void HLTMETAveragePhivsEtaMaps(TFile *f, TString &run_dir)