CMS 3D CMS Logo

JetMETDQMPostProcessor.cc
Go to the documentation of this file.
1 // Migrated to use DQMEDHarvester by: Jyothsna Rani Komaragiri, Oct 2014
2 
4 
9 
10 #include <cmath>
11 #include <cstring>
12 #include <fstream>
13 #include <iomanip>
14 #include <iostream>
15 
17  subDir_ = pset.getUntrackedParameter<std::string>("subDir");
18  patternJetTrg_ = pset.getUntrackedParameter<std::string>("PatternJetTrg", "");
19  patternMetTrg_ = pset.getUntrackedParameter<std::string>("PatternMetTrg", "");
20 }
21 
24  // setup DQM stor //
26 
27  bool isJetDir = false;
28  bool isMetDir = false;
29 
30  TPRegexp patternJet(patternJetTrg_);
31  TPRegexp patternMet(patternMetTrg_);
32 
33  // go to the directory to be processed
34  if (igetter.dirExists(subDir_))
35  ibooker.cd(subDir_);
36  else {
37  edm::LogWarning("JetMETDQMPostProcessor") << "cannot find directory: " << subDir_ << " , skipping";
38  return;
39  }
40 
41  std::vector<std::string> subdirectories = igetter.getSubdirs();
42  for (std::vector<std::string>::iterator dir = subdirectories.begin(); dir != subdirectories.end(); dir++) {
43  ibooker.cd(*dir);
44 
45  isJetDir = false;
46  isMetDir = false;
47 
48  if (TString(*dir).Contains(patternJet))
49  isJetDir = true;
50  if (TString(*dir).Contains(patternMet))
51  isMetDir = true;
52 
53  if (isMetDir) {
54  // std::cout << "JetMETDQMPostProcessor - Met paths: " << ibooker.pwd() <<
55  // " " << *dir << std::endl;
56 
57  // GenMET
58  dividehistos(ibooker,
59  igetter,
60  "_meGenMETTrgMC",
61  "_meGenMET",
62  "_meTurnOngMET",
63  "Gen Missing ET",
64  "Gen Missing ET Turn-On RelVal");
65  dividehistos(ibooker,
66  igetter,
67  "_meGenMETTrg",
68  "_meGenMETTrgLow",
69  "_meTurnOngMETLow",
70  "Gen Missing ETLow",
71  "Gen Missing ET Turn-On Data");
72 
73  // HLTMET
74  dividehistos(ibooker,
75  igetter,
76  "_meHLTMETTrgMC",
77  "_meHLTMET",
78  "_meTurnOnhMET",
79  "HLT Missing ET",
80  "HLT Missing ET Turn-On RelVal");
81  dividehistos(ibooker,
82  igetter,
83  "_meHLTMETTrg",
84  "_meHLTMETTrgLow",
85  "_meTurnOnhMETLow",
86  "HLT Missing ETLow",
87  "HLT Missing ET Turn-On Data");
88  }
89 
90  if (isJetDir) {
91  // std::cout << "JetMETDQMPostProcessor - Jet paths: " << ibooker.pwd() <<
92  // " " << *dir << std::endl;
93 
94  // GenJets
95  dividehistos(ibooker,
96  igetter,
97  "_meGenJetPtTrgMC",
98  "_meGenJetPt",
99  "_meTurnOngJetPt",
100  "Gen Jet Pt",
101  "Gen Jet Pt Turn-On RelVal");
102  dividehistos(ibooker,
103  igetter,
104  "_meGenJetPtTrg",
105  "_meGenJetPtTrgLow",
106  "_meTurnOngJetPt",
107  "Gen Jet PtLow",
108  "Gen Jet Pt Turn-On Data");
109  dividehistos(ibooker,
110  igetter,
111  "_meGenJetEtaTrgMC",
112  "_meGenJetEta",
113  "_meTurnOngJetEta",
114  "Gen Jet Eta",
115  "Gen Jet Eta Turn-On RelVal");
116  dividehistos(ibooker,
117  igetter,
118  "_meGenJetEtaTrg",
119  "_meGenJetEtaTrgLow",
120  "_meTurnOngJetEta",
121  "Gen Jet EtaLow",
122  "Gen Jet Eta Turn-On Data");
123  dividehistos(ibooker,
124  igetter,
125  "_meGenJetPhiTrgMC",
126  "_meGenJetPhi",
127  "_meTurnOngJetPhi",
128  "Gen Jet Phi",
129  "Gen Jet Phi Turn-On RelVal");
130  dividehistos(ibooker,
131  igetter,
132  "_meGenJetPhiTrg",
133  "_meGenJetPhiTrgLow",
134  "_meTurnOngJetPhi",
135  "Gen Jet PhiLow",
136  "Gen Jet Phi Turn-On Data");
137 
138  // HLTJets
139  dividehistos(ibooker,
140  igetter,
141  "_meHLTJetPtTrgMC",
142  "_meHLTJetPt",
143  "_meTurnOnhJetPt",
144  "HLT Jet Pt",
145  "HLT Jet Pt Turn-On RelVal");
146  dividehistos(ibooker,
147  igetter,
148  "_meHLTJetPtTrg",
149  "_meHLTJetPtTrgLow",
150  "_meTurnOnhJetPt",
151  "HLT Jet PtLow",
152  "HLT Jet Pt Turn-On Data");
153  dividehistos(ibooker,
154  igetter,
155  "_meHLTJetEtaTrgMC",
156  "_meHLTJetEta",
157  "_meTurnOnhJetEta",
158  "HLT Jet Eta",
159  "HLT Jet Eta Turn-On RelVal");
160  dividehistos(ibooker,
161  igetter,
162  "_meHLTJetEtaTrg",
163  "_meHLTJetEtaTrgLow",
164  "_meTurnOnhJetEta",
165  "HLT Jet EtaLow",
166  "HLT Jet Eta Turn-On Data");
167  dividehistos(ibooker,
168  igetter,
169  "_meHLTJetPhiTrgMC",
170  "_meHLTJetPhi",
171  "_meTurnOnhJetPhi",
172  "HLT Jet Phi",
173  "HLT Jet Phi Turn-On RelVal");
174  dividehistos(ibooker,
175  igetter,
176  "_meHLTJetPhiTrg",
177  "_meHLTJetPhiTrgLow",
178  "_meTurnOnhJetPhi",
179  "HLT Jet PhiLow",
180  "HLT Jet Phi Turn-On Data");
181  }
182 
183  ibooker.goUp();
184  }
185 }
186 
187 //----------------------------------------------------------------------
189  DQMStore::IGetter &igetter,
190  const std::string &numName,
191  const std::string &denomName,
192  const std::string &outName,
193  const std::string &label,
194  const std::string &titel) {
195  // ibooker.pwd();
196  // std::cout << "In dividehistos: " << ibooker.pwd() << std::endl;
197 
198  // std::cout << numName <<std::endl;
199  TH1F *num = getHistogram(ibooker, igetter, ibooker.pwd() + "/" + numName);
200 
201  // std::cout << denomName << std::endl;
202  TH1F *denom = getHistogram(ibooker, igetter, ibooker.pwd() + "/" + denomName);
203 
204  if (num == nullptr)
205  edm::LogWarning("JetMETDQMPostProcessor")
206  << "numerator histogram " << ibooker.pwd() + "/" + numName << " does not exist";
207  if (denom == nullptr)
208  edm::LogWarning("JetMETDQMPostProcessor")
209  << "denominator histogram " << ibooker.pwd() + "/" + denomName << " does not exist";
210 
211  // Check if histograms actually exist
212  if (!num || !denom)
213  return nullptr;
214 
215  MonitorElement *meOut = ibooker.bookProfile(
216  outName, titel, num->GetXaxis()->GetNbins(), num->GetXaxis()->GetXmin(), num->GetXaxis()->GetXmax(), 0., 1.2);
217  meOut->setEfficiencyFlag();
218  TProfile *out = meOut->getTProfile();
219  out->GetXaxis()->SetTitle(label.c_str());
220  out->SetYTitle("Efficiency");
221  out->SetOption("PE");
222  out->SetLineColor(2);
223  out->SetLineWidth(2);
224  out->SetMarkerStyle(20);
225  out->SetMarkerSize(0.8);
226  out->SetStats(kFALSE);
227 
228  for (int i = 1; i <= num->GetNbinsX(); i++) {
229  double e, low, high;
230  Efficiency((int)num->GetBinContent(i), (int)denom->GetBinContent(i), 0.683, e, low, high);
231  double err = e - low > high - e ? e - low : high - e;
232  // here is the trick to store info in TProfile:
233  out->SetBinContent(i, e);
234  out->SetBinEntries(i, 1);
235  out->SetBinError(i, sqrt(e * e + err * err));
236  }
237  return out;
238 }
239 
240 //----------------------------------------------------------------------
242  DQMStore::IGetter &igetter,
243  const std::string &histoPath) {
244  ibooker.pwd();
245  MonitorElement *monElement = igetter.get(histoPath);
246  if (monElement != nullptr)
247  return monElement->getTH1F();
248  else
249  return nullptr;
250 }
251 
252 //----------------------------------------------------------------------
254  int passing, int total, double level, double &mode, double &lowerBound, double &upperBound) {
255  // protection
256  if (total == 0) {
257  mode = 0.5;
258  lowerBound = 0;
259  upperBound = 1;
260  return;
261  }
262  mode = passing / ((double)total);
263 
264  // see http://root.cern.ch/root/html/TEfficiency.html#compare
265  lowerBound = TEfficiency::Wilson(total, passing, level, false);
266  upperBound = TEfficiency::Wilson(total, passing, level, true);
267 }
268 
269 //----------------------------------------------------------------------
270 
TProfile * getTProfile() const
T getUntrackedParameter(std::string const &, T const &) const
void Efficiency(int passing, int total, double level, double &mode, double &lowerBound, double &upperBound)
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:113
TH1F * getTH1F() const
JetMETDQMPostProcessor(const edm::ParameterSet &pset)
char const * label
void setEfficiencyFlag()
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T sqrt(T t)
Definition: SSEVec.h:18
TH1F * getHistogram(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter, const std::string &histoPath)
TProfile * dividehistos(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter, const std::string &numName, const std::string &denomName, const std::string &outName, const std::string &label, const std::string &titel)
MonitorElement * get(std::string const &path)
Definition: DQMStore.cc:303
std::string const & pwd()
Definition: DQMStore.cc:278
bool dirExists(std::string const &path)
Definition: DQMStore.cc:343
std::vector< std::string > getSubdirs()
Definition: DQMStore.cc:325
dbl *** dir
Definition: mlp_gen.cc:35
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override