CMS 3D CMS Logo

JetMETHLTOfflineClient.cc
Go to the documentation of this file.
1 // Migrated to use DQMEDHarvester by: Jyothsna Rani Komaragiri, Oct 2014
2 
4 
9 
12 
14 {
15  debug_ = false;
16  verbose_ = false;
17 
18  processname_ = iConfig.getParameter<std::string>("processname");
19 
20  hltTag_ = iConfig.getParameter<std::string>("hltTag");
21  if (debug_) std::cout << hltTag_ << std::endl;
22 
23  dirName_=iConfig.getParameter<std::string>("DQMDirName");
24 
25 }
26 
28 
30 {
31  ibooker.setCurrentFolder(dirName_);
32 
33  LogDebug("JetMETHLTOfflineClient") << "dqmEndJob" << std::endl;
34  if (debug_) std::cout << "dqmEndJob" << std::endl;
35 
36  std::vector<MonitorElement*> hltMEs;
37 
38  // Look at all folders, go to the subfolder which includes the string "Eff"
39  std::vector<std::string> fullPathHLTFolders = igetter.getSubdirs();
40  for(auto & fullPathHLTFolder : fullPathHLTFolders) {
41 
42  // Move on only if the folder name contains "Eff" Or "Trigger Summary"
43  if (debug_) std::cout << fullPathHLTFolder << std::endl;
44  if ((fullPathHLTFolder.find("Eff")!=std::string::npos)) {
45  ibooker.setCurrentFolder(fullPathHLTFolder);
46  }
47  else {
48  continue;
49  }
50 
51  // Look at all subfolders, go to the subfolder which includes the string "Eff"
52  std::vector<std::string> fullSubPathHLTFolders = igetter.getSubdirs();
53  for(auto & fullSubPathHLTFolder : fullSubPathHLTFolders) {
54 
55  if (debug_) std::cout << fullSubPathHLTFolder << std::endl;
56  ibooker.setCurrentFolder(fullSubPathHLTFolder);
57 
58  // Look at all MonitorElements in this folder
59  hltMEs = igetter.getContents(fullSubPathHLTFolder);
60  LogDebug("JetMETHLTOfflineClient")<< "Number of MEs for this HLT path = " << hltMEs.size() << std::endl;
61 
62  for(unsigned int k=0;k<hltMEs.size();k++) {
63  if (debug_) std::cout << hltMEs[k]->getName() << std::endl;
64 
65  //-----
66  if ((hltMEs[k]->getName().find("ME_Numerator")!=std::string::npos) && hltMEs[k]->getName().find("ME_Numerator")==0){
67 
68  std::string name = hltMEs[k]->getName();
69  name.erase(0,12); // Removed "ME_Numerator"
70  if (debug_) std::cout <<"==name=="<< name << std::endl;
71  // if( name.find("EtaPhi") !=std::string::npos ) continue; // do not consider EtaPhi 2D plots
72 
73  for(unsigned int l=0;l<hltMEs.size();l++) {
74  if (hltMEs[l]->getName() == "ME_Denominator"+name){
75  // found denominator too
76  if(name.find("EtaPhi") !=std::string::npos) {
77  TH2F* tNumerator = hltMEs[k]->getTH2F();
78  TH2F* tDenominator = hltMEs[l]->getTH2F();
79 
80  std::string title = "Eff_"+hltMEs[k]->getTitle();
81 
82  auto *teff = (TH2F*) tNumerator->Clone(title.c_str());
83  teff->Divide(tNumerator,tDenominator,1,1);
84  ibooker.book2D("ME_Eff_"+name,teff);
85  delete teff;
86  }
87  else{
88  TH1F* tNumerator = hltMEs[k]->getTH1F();
89  TH1F* tDenominator = hltMEs[l]->getTH1F();
90 
91  std::string title = "Eff_"+hltMEs[k]->getTitle();
92 
93  auto *teff = (TH1F*) tNumerator->Clone(title.c_str());
94  teff->Divide(tNumerator,tDenominator,1,1);
95  ibooker.book1D("ME_Eff_"+name,teff);
96  delete teff;
97  }
98  } // Denominator
99  } // Loop-l
100  } // Numerator
101 
102 
103  } // Loop-k
104  } // fullSubPath
105  } // fullPath
106 }
107 
#define LogDebug(id)
std::vector< MonitorElement * > getContents(Args &&...args)
Definition: DQMStore.h:305
T getParameter(std::string const &) const
TH1F * getTH1F() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
int k[5][pyjets_maxn]
~JetMETHLTOfflineClient() override
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:279
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
std::vector< std::string > getSubdirs()
Definition: DQMStore.cc:325
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
JetMETHLTOfflineClient(const edm::ParameterSet &)