CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTJetMETValidation.cc
Go to the documentation of this file.
1 // Migrated to use DQMEDAnalyzer by: Jyothsna Rani Komaragiri, Oct 2014
2 
4 #include "Math/GenVector/VectorUtil.h"
6 
8  triggerEventObject_(consumes<TriggerEventWithRefs>(ps.getUntrackedParameter<edm::InputTag>("triggerEventObject"))),
9  PFJetAlgorithm( consumes<PFJetCollection>(ps.getUntrackedParameter<edm::InputTag>( "PFJetAlgorithm" ) )),
10  GenJetAlgorithm( consumes<GenJetCollection>(ps.getUntrackedParameter<edm::InputTag>( "GenJetAlgorithm" ) )),
11  CaloMETColl( consumes<CaloMETCollection>(ps.getUntrackedParameter<edm::InputTag>( "CaloMETCollection" ) )),
12  GenMETColl( consumes<GenMETCollection>(ps.getUntrackedParameter<edm::InputTag>( "GenMETCollection" ) )),
13  HLTriggerResults(consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>( "HLTriggerResults" ) )),
14  triggerTag_(ps.getUntrackedParameter<std::string>("DQMFolder","SingleJet")),
15  patternJetTrg_(ps.getUntrackedParameter<std::string>("PatternJetTrg","")),
16  patternMetTrg_(ps.getUntrackedParameter<std::string>("PatternMetTrg","")),
17  patternMuTrg_(ps.getUntrackedParameter<std::string>("PatternMuTrg","")),
18  HLTinit_(false)
19 {
20  evtCnt=0;
21 }
22 
24 {
25 }
26 
27 //
28 // member functions
29 //
30 
31 // ------------ method called when starting to processes a run ------------
32 void
34 
35  bool foundMuTrg = false;
36  std::string trgMuNm;
37  bool changedConfig = true;
38 
39  //--define search patterns
40  TPRegexp patternJet(patternJetTrg_);
41  TPRegexp patternMet(patternMetTrg_);
42  TPRegexp patternMu(patternMuTrg_);
43 
44  if (!hltConfig_.init(iRun, iSetup, "HLT", changedConfig)) {
45  edm::LogError("HLTJetMETValidation") << "Initialization of HLTConfigProvider failed!!";
46  return;
47  }
48 
49  std::vector<std::string> validTriggerNames = hltConfig_.triggerNames();
50  for (size_t j = 0; j < validTriggerNames.size(); j++) {
51  //---find the muon path
52  if (TString(validTriggerNames[j]).Contains(patternMu)) {
53  //std::cout <<validTriggerNames[j].c_str()<<std::endl;
54  if (!foundMuTrg) trgMuNm = validTriggerNames[j].c_str();
55  foundMuTrg = true;
56  }
57  //---find the jet paths
58  if (TString(validTriggerNames[j]).Contains(patternJet)) {
59  hltTrgJet.push_back(validTriggerNames[j]);
60  }
61  //---find the met paths
62  if (TString(validTriggerNames[j]).Contains(patternMet)) {
63  hltTrgMet.push_back(validTriggerNames[j]);
64  }
65  }
66 
67  //----set the denominator paths
68  for (size_t it=0;it<hltTrgJet.size();it++) {
69  if (it==0 && foundMuTrg) hltTrgJetLow.push_back(trgMuNm);//--lowest threshold uses muon path
70  if (it==0 && !foundMuTrg) hltTrgJetLow.push_back(hltTrgJet[it]);//---if no muon then itself
71  if (it!=0) hltTrgJetLow.push_back(hltTrgJet[it-1]);
72  //std::cout<<hltTrgJet[it].c_str()<<" "<<hltTrgJetLow[it].c_str()<<std::endl;
73  }
74  int itm(0), itpm(0), itmh(0), itpmh(0);
75  for (size_t it=0;it<hltTrgMet.size();it++) {
76  if (TString(hltTrgMet[it]).Contains("PF")) {
77  if (TString(hltTrgMet[it]).Contains("MHT")) {
78  if( 0 == itpmh ) {
79  if( foundMuTrg ) hltTrgMetLow.push_back(trgMuNm);
80  else hltTrgMetLow.push_back(hltTrgMet[it]);
81  }
82  else hltTrgMetLow.push_back(hltTrgMet[it-1]);
83  itpmh++;
84  }
85  if (TString(hltTrgMet[it]).Contains("MET")) {
86  if( 0 == itpm ) {
87  if( foundMuTrg ) hltTrgMetLow.push_back(trgMuNm);
88  else hltTrgMetLow.push_back(hltTrgMet[it]);
89  }
90  else hltTrgMetLow.push_back(hltTrgMet[it-1]);
91  itpm++;
92  }
93  }
94  else {
95  if (TString(hltTrgMet[it]).Contains("MHT")) {
96  if( 0 == itmh ) {
97  if( foundMuTrg ) hltTrgMetLow.push_back(trgMuNm);
98  else hltTrgMetLow.push_back(hltTrgMet[it]);
99  }
100  else hltTrgMetLow.push_back(hltTrgMet[it-1]);
101  itmh++;
102  }
103  if (TString(hltTrgMet[it]).Contains("MET")) {
104  if( 0 == itm ) {
105  if( foundMuTrg ) hltTrgMetLow.push_back(trgMuNm);
106  else hltTrgMetLow.push_back(hltTrgMet[it]);
107  }
108  else hltTrgMetLow.push_back(hltTrgMet[it-1]);
109  itm++;
110  }
111  }
112  //std::cout<<hltTrgMet[it].c_str()<<" "<<hltTrgMetLow[it].c_str()<<std::endl;
113  }
114 
115 }
116 
117 // ------------ method called to book histograms before starting event loop ------------
118 void
120 {
121 
122  //----define DQM folders and elements
123  for (size_t it=0;it<hltTrgJet.size();it++) {
124  //std::cout<<hltTrgJet[it].c_str()<<" "<<hltTrgJetLow[it].c_str()<<std::endl;
126  //std::cout << "str = " << triggerTag_+hltTrgJet[it].c_str() << std::endl;
127  //std::cout << "trgPathName = " << trgPathName << std::endl;
128  iBooker.setCurrentFolder(trgPathName);
129  _meHLTJetPt.push_back(iBooker.book1D("_meHLTJetPt","Single HLT Jet Pt",100,0,500));
130  _meHLTJetPtTrgMC.push_back(iBooker.book1D("_meHLTJetPtTrgMC","Single HLT Jet Pt - HLT Triggered",100,0,500));
131  _meHLTJetPtTrg.push_back(iBooker.book1D("_meHLTJetPtTrg","Single HLT Jet Pt - HLT Triggered",100,0,500));
132  _meHLTJetPtTrgLow.push_back(iBooker.book1D("_meHLTJetPtTrgLow","Single HLT Jet Pt - HLT Triggered Low",100,0,500));
133 
134  _meHLTJetEta.push_back(iBooker.book1D("_meHLTJetEta","Single HLT Jet Eta",100,-10,10));
135  _meHLTJetEtaTrgMC.push_back(iBooker.book1D("_meHLTJetEtaTrgMC","Single HLT Jet Eta - HLT Triggered",100,-10,10));
136  _meHLTJetEtaTrg.push_back(iBooker.book1D("_meHLTJetEtaTrg","Single HLT Jet Eta - HLT Triggered",100,-10,10));
137  _meHLTJetEtaTrgLow.push_back(iBooker.book1D("_meHLTJetEtaTrgLow","Single HLT Jet Eta - HLT Triggered Low",100,-10,10));
138 
139  _meHLTJetPhi.push_back(iBooker.book1D("_meHLTJetPhi","Single HLT Jet Phi",100,-4.,4.));
140  _meHLTJetPhiTrgMC.push_back(iBooker.book1D("_meHLTJetPhiTrgMC","Single HLT Jet Phi - HLT Triggered",100,-4.,4.));
141  _meHLTJetPhiTrg.push_back(iBooker.book1D("_meHLTJetPhiTrg","Single HLT Jet Phi - HLT Triggered",100,-4.,4.));
142  _meHLTJetPhiTrgLow.push_back(iBooker.book1D("_meHLTJetPhiTrgLow","Single HLT Jet Phi - HLT Triggered Low",100,-4.,4.));
143 
144  _meGenJetPt.push_back(iBooker.book1D("_meGenJetPt","Single Generated Jet Pt",100,0,500));
145  _meGenJetPtTrgMC.push_back(iBooker.book1D("_meGenJetPtTrgMC","Single Generated Jet Pt - HLT Triggered",100,0,500));
146  _meGenJetPtTrg.push_back(iBooker.book1D("_meGenJetPtTrg","Single Generated Jet Pt - HLT Triggered",100,0,500));
147  _meGenJetPtTrgLow.push_back(iBooker.book1D("_meGenJetPtTrgLow","Single Generated Jet Pt - HLT Triggered Low",100,0,500));
148 
149  _meGenJetEta.push_back(iBooker.book1D("_meGenJetEta","Single Generated Jet Eta",100,-10,10));
150  _meGenJetEtaTrgMC.push_back(iBooker.book1D("_meGenJetEtaTrgMC","Single Generated Jet Eta - HLT Triggered",100,-10,10));
151  _meGenJetEtaTrg.push_back(iBooker.book1D("_meGenJetEtaTrg","Single Generated Jet Eta - HLT Triggered",100,-10,10));
152  _meGenJetEtaTrgLow.push_back(iBooker.book1D("_meGenJetEtaTrgLow","Single Generated Jet Eta - HLT Triggered Low",100,-10,10));
153 
154  _meGenJetPhi.push_back(iBooker.book1D("_meGenJetPhi","Single Generated Jet Phi",100,-4.,4.));
155  _meGenJetPhiTrgMC.push_back(iBooker.book1D("_meGenJetPhiTrgMC","Single Generated Jet Phi - HLT Triggered",100,-4.,4.));
156  _meGenJetPhiTrg.push_back(iBooker.book1D("_meGenJetPhiTrg","Single Generated Jet Phi - HLT Triggered",100,-4.,4.));
157  _meGenJetPhiTrgLow.push_back(iBooker.book1D("_meGenJetPhiTrgLow","Single Generated Jet Phi - HLT Triggered Low",100,-4.,4.));
158 
159  }
160  for (size_t it=0;it<hltTrgMet.size();it++) {
161  //std::cout<<hltTrgMet[it].c_str()<<" "<<hltTrgMetLow[it].c_str()<<std::endl;
163  iBooker.setCurrentFolder(trgPathName);
164  _meHLTMET.push_back(iBooker.book1D("_meHLTMET","HLT Missing ET",100,0,500));
165  _meHLTMETTrgMC.push_back(iBooker.book1D("_meHLTMETTrgMC","HLT Missing ET - HLT Triggered",100,0,500));
166  _meHLTMETTrg.push_back(iBooker.book1D("_meHLTMETTrg","HLT Missing ET - HLT Triggered",100,0,500));
167  _meHLTMETTrgLow.push_back(iBooker.book1D("_meHLTMETTrgLow","HLT Missing ET - HLT Triggered Low",100,0,500));
168 
169  _meGenMET.push_back(iBooker.book1D("_meGenMET","Generated Missing ET",100,0,500));
170  _meGenMETTrgMC.push_back(iBooker.book1D("_meGenMETTrgMC","Generated Missing ET - HLT Triggered",100,0,500));
171  _meGenMETTrg.push_back(iBooker.book1D("_meGenMETTrg","Generated Missing ET - HLT Triggered",100,0,500));
172  _meGenMETTrgLow.push_back(iBooker.book1D("_meGenMETTrgLow","Generated Missing ET - HLT Triggered Low",100,0,500));
173  }
174 
175 }
176 
177 // ------------ method called for each event ------------
178 void
180 {
181  using namespace std;
182  using namespace edm;
183  using namespace reco;
184  using namespace l1extra;
185  using namespace trigger;
186 
187  evtCnt++;
188 
189  //get The triggerEvent
191  iEvent.getByToken(triggerEventObject_,trigEv);
192 
193  //get TriggerResults object
194  bool gotHLT=true;
195  std::vector<bool> myTrigJ;
196  for (size_t it=0;it<hltTrgJet.size();it++) myTrigJ.push_back(false);
197  std::vector<bool> myTrigJLow;
198  for (size_t it=0;it<hltTrgJetLow.size();it++) myTrigJLow.push_back(false);
199  std::vector<bool> myTrigM;
200  for (size_t it=0;it<hltTrgMet.size();it++) myTrigM.push_back(false);
201  std::vector<bool> myTrigMLow;
202  for (size_t it=0;it<hltTrgMetLow.size();it++) myTrigMLow.push_back(false);
203 
204 
205  Handle<TriggerResults> hltresults;
206  iEvent.getByToken(HLTriggerResults,hltresults);
207  if (! hltresults.isValid() ) {
208  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No HLTRESULTS";
209  gotHLT=false;
210  }
211 
212  if (gotHLT) {
213  const edm::TriggerNames & triggerNames = iEvent.triggerNames(*hltresults);
214  getHLTResults(*hltresults, triggerNames);
215 
216  //---pick-up the jet trigger decisions
217  for (size_t it=0;it<hltTrgJet.size();it++) {
219  if (trig_iter==hltTriggerMap.end()){
220  //std::cout << "Could not find trigger path with name: " << _probefilter.label() << std::endl;
221  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << "Could not find trigger path with name: " << _probefilter.label();
222  }else{
223  myTrigJ[it]=trig_iter->second;
224  }
225  //std::cout<<hltTrgJet[it].c_str()<<" "<<myTrigJ[it]<<std::endl;
226  }
227  for (size_t it=0;it<hltTrgJetLow.size();it++) {
229  if (trig_iter==hltTriggerMap.end()){
230  //std::cout << "Could not find trigger path with name: " << _probefilter.label() << std::endl;
231  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << "Could not find trigger path with name: " << _probefilter.label();
232  }else{
233  myTrigJLow[it]=trig_iter->second;
234  }
235  //std::cout<<hltTrgJetLow[it].c_str()<<" "<<myTrigJLow[it]<<std::endl;
236  }
237  //---pick-up the met trigger decisions
238  for (size_t it=0;it<hltTrgMet.size();it++) {
240  if (trig_iter==hltTriggerMap.end()){
241  //std::cout << "Could not find trigger path with name: " << _probefilter.label() << std::endl;
242  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << "Could not find trigger path with name: " << _probefilter.label();
243  }else{
244  myTrigM[it]=trig_iter->second;
245  }
246  //std::cout<<hltTrgMet[it].c_str()<<" "<<myTrigM[it]<<std::endl;
247  }
248  for (size_t it=0;it<hltTrgMetLow.size();it++) {
250  if (trig_iter==hltTriggerMap.end()){
251  //std::cout << "Could not find trigger path with name: " << _probefilter.label() << std::endl;
252  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << "Could not find trigger path with name: " << _probefilter.label();
253  }else{
254  myTrigMLow[it]=trig_iter->second;
255  }
256  //std::cout<<hltTrgMetLow[it].c_str()<<" "<<myTrigMLow[it]<<std::endl;
257  }
258  }
259 
260  // --- Fill histos for PFJet paths ---
261  // HLT jets namely hltAK4PFJets
263  iEvent.getByToken( PFJetAlgorithm, pfJets );
264  double pfJetPt=-1.;
265  double pfJetEta=-999.;
266  double pfJetPhi=-999.;
267 
268  if (pfJets.isValid()) {
269  //Loop over the PFJets and fill some histograms
270  int jetInd = 0;
271  for( PFJetCollection::const_iterator pf = pfJets->begin(); pf != pfJets->end(); ++ pf ) {
272  //std::cout << "PF JET #" << jetInd << std::endl << pf->print() << std::endl;
273  if (jetInd == 0){
274  pfJetPt=pf->pt();
275  pfJetEta=pf->eta();
276  pfJetPhi=pf->phi();
277  for (size_t it=0;it<hltTrgJet.size();it++) {
278  _meHLTJetPt[it]->Fill( pfJetPt );
279  _meHLTJetEta[it]->Fill( pfJetEta );
280  _meHLTJetPhi[it]->Fill( pfJetPhi );
281  if (myTrigJ[it]) _meHLTJetPtTrgMC[it]->Fill( pfJetPt );
282  if (myTrigJ[it]) _meHLTJetEtaTrgMC[it]->Fill( pfJetEta );
283  if (myTrigJ[it]) _meHLTJetPhiTrgMC[it]->Fill( pfJetPhi );
284  if (myTrigJ[it] && myTrigJLow[it]) _meHLTJetPtTrg[it]->Fill( pfJetPt );
285  if (myTrigJ[it] && myTrigJLow[it]) _meHLTJetEtaTrg[it]->Fill( pfJetEta );
286  if (myTrigJ[it] && myTrigJLow[it]) _meHLTJetPhiTrg[it]->Fill( pfJetPhi );
287  if (myTrigJLow[it]) _meHLTJetPtTrgLow[it]->Fill( pfJetPt );
288  if (myTrigJLow[it]) _meHLTJetEtaTrgLow[it]->Fill( pfJetEta );
289  if (myTrigJLow[it]) _meHLTJetPhiTrgLow[it]->Fill( pfJetPhi );
290  }
291  jetInd++;
292  }
293  }//loop over pfjets
294  }
295  else{
296  //std::cout << " -- No PFJets" << std::endl;
297  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No PFJets";
298  }
299 
300  //GenJets
301  Handle<GenJetCollection> genJets;
302  iEvent.getByToken( GenJetAlgorithm, genJets );
303  double genJetPt=-1.;
304  double genJetEta=-999.;
305  double genJetPhi=-999.;
306 
307  if (genJets.isValid()) {
308  //Loop over the GenJets and fill some histograms
309  int jetInd = 0;
310  for( GenJetCollection::const_iterator gen = genJets->begin(); gen != genJets->end(); ++ gen ) {
311  if (jetInd == 0){
312  genJetPt=gen->pt();
313  genJetEta=gen->eta();
314  genJetPhi=gen->phi();
315  for (size_t it=0;it<hltTrgJet.size();it++) {
316  _meGenJetPt[it]->Fill( genJetPt );
317  _meGenJetEta[it]->Fill( genJetEta );
318  _meGenJetPhi[it]->Fill( genJetPhi );
319  if (myTrigJ[it]) _meGenJetPtTrgMC[it]->Fill( genJetPt );
320  if (myTrigJ[it]) _meGenJetEtaTrgMC[it]->Fill( genJetEta );
321  if (myTrigJ[it]) _meGenJetPhiTrgMC[it]->Fill( genJetPhi );
322  if (myTrigJ[it] && myTrigJLow[it]) _meGenJetPtTrg[it]->Fill( genJetPt );
323  if (myTrigJ[it] && myTrigJLow[it]) _meGenJetEtaTrg[it]->Fill( genJetEta );
324  if (myTrigJ[it] && myTrigJLow[it]) _meGenJetPhiTrg[it]->Fill( genJetPhi );
325  if (myTrigJLow[it]) _meGenJetPtTrgLow[it]->Fill( genJetPt );
326  if (myTrigJLow[it]) _meGenJetEtaTrgLow[it]->Fill( genJetEta );
327  if (myTrigJLow[it]) _meGenJetPhiTrgLow[it]->Fill( genJetPhi );
328  }
329  jetInd++;
330  }
331  }
332  }
333  else{
334  //std::cout << " -- No GenJets" << std::endl;
335  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No GenJets";
336  }
337 
338  // --- Fill histos for PFMET paths ---
339  // HLT MET namely hltmet
341  iEvent.getByToken(CaloMETColl, recmet);
342 
343  double calMet=-1;
344  if (recmet.isValid()) {
345  typedef CaloMETCollection::const_iterator cmiter;
346  //std::cout << "Size of MET collection" << recmet.size() << std::endl;
347  for ( cmiter i=recmet->begin(); i!=recmet->end(); i++) {
348  calMet = i->pt();
349  for (size_t it=0;it<hltTrgMet.size();it++) {
350  _meHLTMET[it] -> Fill(calMet);
351  if (myTrigM.size() > it && myTrigM[it]) _meHLTMETTrgMC[it] -> Fill(calMet);
352  if (myTrigM.size() > it && myTrigMLow.size() > it && myTrigM[it] && myTrigMLow[it]) _meHLTMETTrg[it] -> Fill(calMet);
353  if (myTrigMLow.size() > it && myTrigMLow[it]) _meHLTMETTrgLow[it] -> Fill(calMet);
354  }
355  }
356  }
357  else{
358  //std::cout << " -- No MET Collection with name: " << CaloMETColl << std::endl;
359  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No MET Collection with name: "<< CaloMETColl;
360  }
361 
363  iEvent.getByToken(GenMETColl, genmet);
364 
365  double genMet=-1;
366  if (genmet.isValid()) {
367  typedef GenMETCollection::const_iterator cmiter;
368  for ( cmiter i=genmet->begin(); i!=genmet->end(); i++) {
369  genMet = i->pt();
370  for (size_t it=0;it<hltTrgMet.size();it++) {
371  _meGenMET[it] -> Fill(genMet);
372  if (myTrigM.size() > it && myTrigM[it]) _meGenMETTrgMC[it] -> Fill(genMet);
373  if (myTrigM.size() > it && myTrigMLow.size() > it && myTrigM[it] && myTrigMLow[it]) _meGenMETTrg[it] -> Fill(genMet);
374  if (myTrigMLow.size() > it && myTrigMLow[it]) _meGenMETTrgLow[it] -> Fill(genMet);
375  }
376  }
377  }
378  else{
379  //std::cout << " -- No GenMET Collection with name: " << GenMETColl << std::endl;
380  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No GenMET Collection with name: "<< GenMETColl;
381  }
382 
383 }
384 
386  const edm::TriggerNames & triggerNames) {
387 
388  int ntrigs=hltresults.size();
389  if (! HLTinit_){
390  HLTinit_=true;
391 
392  for (int itrig = 0; itrig != ntrigs; ++itrig){
393  std::string trigName = triggerNames.triggerName(itrig);
394  // std::cout << "trigger " << itrig << ": " << trigName << std::endl;
395  }
396  }
397 
398  for (int itrig = 0; itrig != ntrigs; ++itrig){
399  std::string trigName = triggerNames.triggerName(itrig);
400  bool accept=hltresults.accept(itrig);
401 
402  //if (accept) _triggerResults->Fill(float(itrig));
403 
404  // fill the trigger map
406  trig_iter=hltTriggerMap.find(trigName);
407  if (trig_iter==hltTriggerMap.end())
408  hltTriggerMap.insert(valType(trigName,accept));
409  else
410  trig_iter->second=accept;
411  }
412 }
std::vector< MonitorElement * > _meHLTMETTrgMC
std::vector< std::string > hltTrgMet
std::vector< MonitorElement * > _meGenJetEtaTrgMC
std::vector< MonitorElement * > _meGenJetPt
int i
Definition: DBlmapReader.cc:9
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:208
Collection of Gen MET.
std::vector< MonitorElement * > _meGenJetPhiTrg
std::vector< std::string > hltTrgJet
std::vector< MonitorElement * > _meHLTJetPtTrg
std::map< std::string, bool > hltTriggerMap
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &iRun, edm::EventSetup const &iSetup) override
edm::EDGetTokenT< TriggerEventWithRefs > triggerEventObject_
InputTag of TriggerEventWithRefs to analyze.
edm::EDGetTokenT< GenMETCollection > GenMETColl
std::vector< MonitorElement * > _meHLTJetPhiTrgLow
bool accept() const
Has at least one path accepted the event?
std::vector< GenJet > GenJetCollection
collection of GenJet objects
const std::vector< std::string > & triggerNames() const
names of trigger paths
edm::EDGetTokenT< GenJetCollection > GenJetAlgorithm
std::vector< MonitorElement * > _meGenMET
virtual void dqmBeginRun(edm::Run const &iRun, edm::EventSetup const &iSetup) override
std::vector< MonitorElement * > _meHLTMETTrgLow
std::vector< MonitorElement * > _meGenJetPtTrgMC
std::vector< MonitorElement * > _meGenMETTrgMC
edm::EDGetTokenT< CaloMETCollection > CaloMETColl
std::vector< std::string > hltTrgMetLow
std::vector< MonitorElement * > _meGenJetPhiTrgMC
std::vector< MonitorElement * > _meGenMETTrgLow
std::vector< MonitorElement * > _meHLTJetPtTrgLow
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:25
HLTJetMETValidation(const edm::ParameterSet &)
std::vector< MonitorElement * > _meHLTJetEta
static const std::string removeVersion(const std::string &trigger)
Collection of Calo MET.
int iEvent
Definition: GenABIO.cc:230
std::vector< MonitorElement * > _meHLTJetPhiTrgMC
edm::EDGetTokenT< edm::TriggerResults > HLTriggerResults
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
unsigned int size() const
Get number of paths stored.
std::vector< MonitorElement * > _meHLTJetEtaTrg
std::vector< MonitorElement * > _meHLTJetEtaTrgLow
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
edm::EDGetTokenT< PFJetCollection > PFJetAlgorithm
int j
Definition: DBlmapReader.cc:9
std::vector< MonitorElement * > _meGenJetEtaTrg
std::vector< MonitorElement * > _meGenJetPtTrg
Container::value_type value_type
std::vector< MonitorElement * > _meGenMETTrg
std::vector< MonitorElement * > _meGenJetEtaTrgLow
std::vector< MonitorElement * > _meHLTJetPt
void getHLTResults(const edm::TriggerResults &, const edm::TriggerNames &triggerNames)
std::vector< MonitorElement * > _meHLTJetPhi
std::vector< MonitorElement * > _meGenJetEta
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:27
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
std::vector< MonitorElement * > _meGenJetPhi
std::vector< MonitorElement * > _meGenJetPtTrgLow
TLorentzVector genMet(const HepMC::GenEvent *all, double etamin=-9999., double etamax=9999.)
HLTConfigProvider hltConfig_
std::vector< MonitorElement * > _meHLTJetPtTrgMC
std::vector< MonitorElement * > _meHLTMETTrg
std::vector< MonitorElement * > _meHLTJetPhiTrg
std::vector< PFJet > PFJetCollection
collection of PFJet objects
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
std::vector< std::string > hltTrgJetLow
volatile std::atomic< bool > shutdown_flag false
std::vector< MonitorElement * > _meGenJetPhiTrgLow
std::map< std::string, bool >::iterator trig_iter
tuple pfJets
Definition: pfJets_cff.py:8
std::vector< MonitorElement * > _meHLTJetEtaTrgMC
Definition: Run.h:41
std::vector< MonitorElement * > _meHLTMET