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