ATLAS Offline Software
Loading...
Searching...
No Matches
MonitoringFile_HLTCaloAveragePtPhiEtaMaps.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 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, std::string& 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.find("run") == std::string::npos) {
36 continue;
37 }
38
39
40
41 // begin HLTMET
42 // note 1: prefix all dirs and hists with '/'
43 // note 2: missing dir => return
44 // note 3: missing hist => continue
45 std::string hlt_top = run_dir + "/HLT"; // toplevel
46
47 std::vector<std::string> calo_fold;
48 std::vector<std::pair<std::string, std::string> > hist_numr;
49
50 // The folders that we want to get the hists from
51 calo_fold.push_back("/CaloMon");
52 calo_fold.push_back("/CaloMonL2");
53
54 // check if the folders are in hlt
55 for (std::vector<std::string>::iterator it = calo_fold.begin(); it != calo_fold.end(); ++it) {
56 std::string theHistDir = hlt_top + *it;
57 TDirectory* dir = f->GetDirectory(theHistDir.c_str());
58 if (!dir) {
59 if (debugLevel > 0) std::cerr << "--> HLTCaloAveragePtPhiEtaMaps: directory " << theHistDir << " not found" << std::endl;
60 return;
61 }
62 }
63
64 // pairs of Num and Dem
65 hist_numr.push_back(std::make_pair<std::string, std::string>("/EnergyAccetaphiLAr", "/HitAccetaphiLAr"));
66 hist_numr.push_back(std::make_pair<std::string, std::string>("/EnergyAccetaphiTile", "/HitAccetaphiTile"));
67
68
69 // we have all dirs, get the component histograms
70 for (std::vector<std::string>::iterator itFex = calo_fold.begin(); itFex != calo_fold.end(); ++itFex) {
71 for (std::vector<std::pair<std::string, std::string> >::iterator itNum = hist_numr.begin(); itNum != hist_numr.end(); ++itNum) {
72 TH2F* hnum(0), *hden(0);
73
74 // prepend histogram name with path and append with suffix [_00 .., _24 for each component]
75 std::string thePath = hlt_top + (*itFex);
76 std::string numPath = thePath + (*itNum).first;
77 std::string denPath = thePath + (*itNum).second;
78
79 // get histograms
80 hnum = (TH2F*) (f->Get(numPath.c_str()));
81 hden = (TH2F*) (f->Get(denPath.c_str()));
82 // test if histograms are present
83 if (!hnum) {
84 if (debugLevel > 0) std::cerr << "--> HLTCaloAveragePtPhiEtaMaps: no histogram " << numPath << std::endl;
85 continue;
86 }
87 if (!hden) {
88 if (debugLevel > 0) std::cerr << "--> HLTCaloAveragePtPhiEtaMaps: no histogram " << denPath << std::endl;
89 continue;
90 }
91
92
93 // get directory of histograms
94 TDirectory* dir = f->GetDirectory(thePath.c_str());
95
96 // these are disabled, because we have to worry about updating metadata
97 // clone numerator histogram in the same directory; prepend with "avg_"
98 // use numerator to do the job
99
100 // divide num by den to get average quantity.
101 hnum->Divide(hnum, hden);
102
103 // fix histogram titles
104 std::string title = hnum->GetTitle();
105 hnum->SetTitle("Average Transverse Energy per eta/phi bin");
106
107 dir->cd();
108 hnum->Write("", TObject::kOverwrite);
109 } // done looping over quantities
110 } // done looping over directories
111 } // end while loop over all keys
112 } // end method MonitoringFile::HLTCaloAveragePtPhiEtaMaps
113} // end namespace
static int getDebugLevel()
static void HLTCaloAveragePtPhiEtaMaps(TFile *f, std::string &run_dir)