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.
2 #include "Math/GenVector/VectorUtil.h"
4 
6  triggerEventObject_(ps.getUntrackedParameter<edm::InputTag>("triggerEventObject")),
7  CaloJetAlgorithm( ps.getUntrackedParameter<edm::InputTag>( "CaloJetAlgorithm" ) ),
8  GenJetAlgorithm( ps.getUntrackedParameter<edm::InputTag>( "GenJetAlgorithm" ) ),
9  CaloMETColl( ps.getUntrackedParameter<edm::InputTag>( "CaloMETCollection" ) ),
10  GenMETColl( ps.getUntrackedParameter<edm::InputTag>( "GenMETCollection" ) ),
11  HLTriggerResults( ps.getParameter<edm::InputTag>( "HLTriggerResults" ) ),
12  triggerTag_(ps.getUntrackedParameter<std::string>("DQMFolder","SingleJet")),
13  _HLTPath(ps.getUntrackedParameter<edm::InputTag>("HLTPath")),
14  _HLTLow(ps.getUntrackedParameter<edm::InputTag>("HLTLow")),
15  outFile_(ps.getUntrackedParameter<std::string>("OutputFileName","")),
16  HLTinit_(false),
17  //JL
18  writeFile_(ps.getUntrackedParameter<bool>("WriteFile",false))
19 {
20  evtCnt=0;
21 
22  //Declare DQM Store
23  DQMStore* store = &*edm::Service<DQMStore>();
24 
25  if(store)
26  {
27  //Create the histograms
29  _meRecoJetPt= store->book1D("_meRecoJetPt","Single Reconstructed Jet Pt",100,0,500);
30  _meRecoJetPtTrg= store->book1D("_meRecoJetPtTrg","Single Reconstructed Jet Pt -- HLT Triggered",100,0,500);
31  _meRecoJetPtTrgLow= store->book1D("_meRecoJetPtTrgLow","Single Reconstructed Jet Pt -- HLT Triggered Low",100,0,500);
32 
33  _meRecoJetEta= store->book1D("_meRecoJetEta","Single Reconstructed Jet Eta",100,-10,10);
34  _meRecoJetEtaTrg= store->book1D("_meRecoJetEtaTrg","Single Reconstructed Jet Eta -- HLT Triggered",100,-10,10);
35  _meRecoJetEtaTrgLow= store->book1D("_meRecoJetEtaTrgLow","Single Reconstructed Jet Eta -- HLT Triggered Low",100,-10,10);
36 
37  _meRecoJetPhi= store->book1D("_meRecoJetPhi","Single Reconstructed Jet Phi",100,-4.,4.);
38  _meRecoJetPhiTrg= store->book1D("_meRecoJetPhiTrg","Single Reconstructed Jet Phi -- HLT Triggered",100,-4.,4.);
39  _meRecoJetPhiTrgLow= store->book1D("_meRecoJetPhiTrgLow","Single Reconstructed Jet Phi -- HLT Triggered Low",100,-4.,4.);
40 
41  _meGenJetPt= store->book1D("_meGenJetPt","Single Generated Jet Pt",100,0,500);
42  _meGenJetPtTrg= store->book1D("_meGenJetPtTrg","Single Generated Jet Pt -- HLT Triggered",100,0,500);
43  _meGenJetPtTrgLow= store->book1D("_meGenJetPtTrgLow","Single Generated Jet Pt -- HLT Triggered Low",100,0,500);
44 
45  _meGenJetEta= store->book1D("_meGenJetEta","Single Generated Jet Eta",100,-10,10);
46  _meGenJetEtaTrg= store->book1D("_meGenJetEtaTrg","Single Generated Jet Eta -- HLT Triggered",100,-10,10);
47  _meGenJetEtaTrgLow= store->book1D("_meGenJetEtaTrgLow","Single Generated Jet Eta -- HLT Triggered Low",100,-10,10);
48 
49  _meGenJetPhi= store->book1D("_meGenJetPhi","Single Generated Jet Phi",100,-4.,4.);
50  _meGenJetPhiTrg= store->book1D("_meGenJetPhiTrg","Single Generated Jet Phi -- HLT Triggered",100,-4.,4.);
51  _meGenJetPhiTrgLow= store->book1D("_meGenJetPhiTrgLow","Single Generated Jet Phi -- HLT Triggered Low",100,-4.,4.);
52 
53  _meRecoMET= store->book1D("_meRecoMET","Reconstructed Missing ET",100,0,500);
54  _meRecoMETTrg= store->book1D("_meRecoMETTrg","Reconstructed Missing ET -- HLT Triggered",100,0,500);
55  _meRecoMETTrgLow= store->book1D("_meRecoMETTrgLow","Reconstructed Missing ET -- HLT Triggered Low",100,0,500);
56 
57  _meGenMET= store->book1D("_meGenMET","Generated Missing ET",100,0,500);
58  _meGenMETTrg= store->book1D("_meGenMETTrg","Generated Missing ET -- HLT Triggered",100,0,500);
59  _meGenMETTrgLow= store->book1D("_meGenMETTrgLow","Generated Missing ET -- HLT Triggered Low",100,0,500);
60 
61  _meGenHT= store->book1D("_meGenHT","Generated HT",100,0,1000);
62  _meGenHTTrg= store->book1D("_meGenHTTrg","Generated HT -- HLT Triggered",100,0,1000);
63  _meGenHTTrgLow= store->book1D("_meGenHTTrgLow","Generated HT -- HLT Triggered Low",100,0,1000);
64 
65  _meRecoHT= store->book1D("_meRecoHT","Reconstructed HT",100,0,1000);
66  _meRecoHTTrg= store->book1D("_meRecoHTTrg","Reconstructed HT -- HLT Triggered",100,0,1000);
67  _meRecoHTTrgLow= store->book1D("_meRecoHTTrgLow","Reconstructed HT -- HLT Triggered Low",100,0,1000);
68 
69  _triggerResults = store->book1D( "_triggerResults", "HLT Results", 200, 0, 200 );
70 
71  }
72 
73  //if (writeFile_) printf("Initializing\n");
74 }
75 
77 {
78 }
79 
80 //
81 // member functions
82 //
83 
84 void
86 {
87 
88  //Write DQM thing..
89  if(outFile_.size()>0)
92  }
93 
94 }
95 
96 void
98 {
99  using namespace std;
100  using namespace edm;
101  using namespace reco;
102  using namespace l1extra;
103  using namespace trigger;
104 
105  evtCnt++;
106  //get The triggerEvent
107 
109  iEvent.getByLabel(triggerEventObject_,trigEv);
110 
111 // get TriggerResults object
112 
113  bool gotHLT=true;
114  bool myTrig=false;
115  bool myTrigLow=false;
116 
117  Handle<TriggerResults> hltresults,hltresultsDummy;
118  iEvent.getByLabel(HLTriggerResults,hltresults);
119  if (! hltresults.isValid() ) {
120  //std::cout << " -- No HLTRESULTS";
121  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No HLTRESULTS";
122  gotHLT=false;
123  }
124 
125  if (gotHLT) {
126  const edm::TriggerNames & triggerNames = iEvent.triggerNames(*hltresults);
127  getHLTResults(*hltresults, triggerNames);
128  // trig_iter=hltTriggerMap.find(MyTrigger);
130  if (trig_iter==hltTriggerMap.end()){
131  //std::cout << "Could not find trigger path with name: " << _probefilter.label() << std::endl;
132  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << "Could not find trigger path with name: " << _probefilter.label();
133  }else{
134  myTrig=trig_iter->second;
135  }
137  if (trig_iter==hltTriggerMap.end()){
138  //std::cout << "Could not find trigger path with name: " << _probefilter.label() << std::endl;
139  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << "Could not find trigger path with name: " << _probefilter.label();
140  }else{
141  myTrigLow=trig_iter->second;
142  }
143  }
144 
145  Handle<CaloJetCollection> caloJets,caloJetsDummy;
146  iEvent.getByLabel( CaloJetAlgorithm, caloJets );
147  double calJetPt=-1.;
148  double calJetEta=-999.;
149  double calJetPhi=-999.;
150  double calHT=0;
151  if (caloJets.isValid()) {
152  //Loop over the CaloJets and fill some histograms
153  int jetInd = 0;
154  for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end(); ++ cal ) {
155  // std::cout << "CALO JET #" << jetInd << std::endl << cal->print() << std::endl;
156  // h_ptCal->Fill( cal->pt() );
157  if (jetInd == 0){
158  //h_ptCalLeading->Fill( cal->pt() );
159  calJetPt=cal->pt();
160  calJetEta=cal->eta();
161  calJetPhi=cal->phi();
162  _meRecoJetPt->Fill( calJetPt );
163  _meRecoJetEta->Fill( calJetEta );
164  _meRecoJetPhi->Fill( calJetPhi );
165  if (myTrig) _meRecoJetPtTrg->Fill( calJetPt );
166  if (myTrig) _meRecoJetEtaTrg->Fill( calJetEta );
167  if (myTrig) _meRecoJetPhiTrg->Fill( calJetPhi );
168  if (myTrigLow) _meRecoJetPtTrgLow->Fill( calJetPt );
169  if (myTrigLow) _meRecoJetEtaTrgLow->Fill( calJetEta );
170  if (myTrigLow) _meRecoJetPhiTrgLow->Fill( calJetPhi );
171 
172  //h_etaCalLeading->Fill( cal->eta() );
173  //h_phiCalLeading->Fill( cal->phi() );
174 
175 
176  jetInd++;
177  }
178  if (cal->pt()>30) {
179  calHT+=cal->pt();
180  }
181  }
182  _meRecoHT->Fill( calHT );
183  if (myTrig) _meRecoHTTrg->Fill( calHT );
184  if (myTrigLow) _meRecoHTTrgLow->Fill( calHT );
185  }else{
186  //std::cout << " -- No CaloJets" << std::endl;
187  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No CaloJets";
188  }
189 
190  Handle<GenJetCollection> genJets,genJetsDummy;
191  iEvent.getByLabel( GenJetAlgorithm, genJets );
192  double genJetPt=-1.;
193  double genJetEta=-999.;
194  double genJetPhi=-999.;
195  double genHT=0;
196  if (genJets.isValid()) {
197  //Loop over the GenJets and fill some histograms
198  int jetInd = 0;
199  for( GenJetCollection::const_iterator gen = genJets->begin(); gen != genJets->end(); ++ gen ) {
200  // std::cout << "CALO JET #" << jetInd << std::endl << cal->print() << std::endl;
201  // h_ptCal->Fill( cal->pt() );
202  if (jetInd == 0){
203  //h_ptCalLeading->Fill( gen->pt() );
204  genJetPt=gen->pt();
205  genJetEta=gen->eta();
206  genJetPhi=gen->phi();
207  _meGenJetPt->Fill( genJetPt );
208  _meGenJetEta->Fill( genJetEta );
209  _meGenJetPhi->Fill( genJetPhi );
210  if (myTrig) _meGenJetPtTrg->Fill( genJetPt );
211  if (myTrig) _meGenJetEtaTrg->Fill( genJetEta );
212  if (myTrig) _meGenJetPhiTrg->Fill( genJetPhi );
213  if (myTrigLow) _meGenJetPtTrgLow->Fill( genJetPt );
214  if (myTrigLow) _meGenJetEtaTrgLow->Fill( genJetEta );
215  if (myTrigLow) _meGenJetPhiTrgLow->Fill( genJetPhi );
216  //h_etaCalLeading->Fill( gen->eta() );
217  //h_phiCalLeading->Fill( gen->phi() );
218 
219  //if (myTrig) h_ptCalTrig->Fill( gen->pt() );
220  jetInd++;
221  }
222  if (gen->pt()>30) {
223  genHT+=gen->pt();
224  }
225  }
226  _meGenHT->Fill( genHT );
227  if (myTrig) _meGenHTTrg->Fill( genHT );
228  if (myTrigLow) _meGenHTTrgLow->Fill( genHT );
229  }else{
230  //std::cout << " -- No GenJets" << std::endl;
231  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No GenJets";
232  }
233 
234  edm::Handle<CaloMETCollection> recmet, recmetDummy;
235  iEvent.getByLabel(CaloMETColl,recmet);
236 
237  double calMet=-1;
238  if (recmet.isValid()) {
239  typedef CaloMETCollection::const_iterator cmiter;
240  //std::cout << "Size of MET collection" << recmet.size() << std::endl;
241  for ( cmiter i=recmet->begin(); i!=recmet->end(); i++) {
242  calMet = i->pt();
243  //mcalphi = i->phi();
244  //mcalsum = i->sumEt();
245  _meRecoMET -> Fill(calMet);
246  if (myTrig) _meRecoMETTrg -> Fill(calMet);
247  if (myTrigLow) _meRecoMETTrgLow -> Fill(calMet);
248  }
249  }else{
250  //std::cout << " -- No MET Collection with name: " << CaloMETColl << std::endl;
251  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No MET Collection with name: "<< CaloMETColl;
252  }
253 
254  edm::Handle<GenMETCollection> genmet, genmetDummy;
255  iEvent.getByLabel(GenMETColl,genmet);
256 
257  double genMet=-1;
258  if (genmet.isValid()) {
259  typedef GenMETCollection::const_iterator cmiter;
260  //std::cout << "Size of GenMET collection" << recmet.size() << std::endl;
261  for ( cmiter i=genmet->begin(); i!=genmet->end(); i++) {
262  genMet = i->pt();
263  //mcalphi = i->phi();
264  //mcalsum = i->sumEt();
265  _meGenMET -> Fill(genMet);
266  if (myTrig) _meGenMETTrg -> Fill(genMet);
267  if (myTrigLow) _meGenMETTrgLow -> Fill(genMet);
268  }
269  }else{
270  //std::cout << " -- No GenMET Collection with name: " << GenMETColl << std::endl;
271  //if (evtCnt==1) edm::LogWarning("HLTJetMETValidation") << " -- No GenMET Collection with name: "<< GenMETColl;
272  }
273 }
274 
276  const edm::TriggerNames & triggerNames) {
277 
278  int ntrigs=hltresults.size();
279  if (! HLTinit_){
280  HLTinit_=true;
281 
282  //if (writeFile_) std::cout << "Number of HLT Paths: " << ntrigs << std::endl;
283  for (int itrig = 0; itrig != ntrigs; ++itrig){
284  std::string trigName = triggerNames.triggerName(itrig);
285  // std::cout << "trigger " << itrig << ": " << trigName << std::endl;
286  }
287  }
288 
289  for (int itrig = 0; itrig != ntrigs; ++itrig){
290  std::string trigName = triggerNames.triggerName(itrig);
291  bool accept=hltresults.accept(itrig);
292 
293  if (accept) _triggerResults->Fill(float(itrig));
294 
295  // fill the trigger map
297  trig_iter=hltTriggerMap.find(trigName);
298  if (trig_iter==hltTriggerMap.end())
299  hltTriggerMap.insert(valType(trigName,accept));
300  else
301  trig_iter->second=accept;
302  }
303 }
MonitorElement * _meGenJetEtaTrgLow
MonitorElement * _meRecoJetPtTrgLow
int i
Definition: DBlmapReader.cc:9
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:208
MonitorElement * _meGenJetEta
MonitorElement * _meRecoJetEtaTrgLow
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:519
MonitorElement * _meRecoJetPtTrg
std::map< std::string, bool > hltTriggerMap
MonitorElement * _triggerResults
MonitorElement * _meRecoJetPhiTrgLow
MonitorElement * _meRecoJetEta
MonitorElement * _meRecoJetEtaTrg
MonitorElement * _meGenJetPt
bool accept() const
Has at least one path accepted the event?
MonitorElement * _meRecoJetPhiTrg
MonitorElement * _meRecoHTTrgLow
MonitorElement * _meGenMETTrgLow
MonitorElement * _meRecoMETTrg
Container::value_type value_type
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:21
HLTJetMETValidation(const edm::ParameterSet &)
MonitorElement * _meGenJetPtTrgLow
MonitorElement * _meRecoMET
void Fill(long long x)
edm::InputTag HLTriggerResults
int iEvent
Definition: GenABIO.cc:243
MonitorElement * _meRecoHTTrg
unsigned int size() const
Get number of paths stored.
edm::InputTag GenJetAlgorithm
MonitorElement * _meGenHTTrgLow
MonitorElement * _meGenHT
MonitorElement * _meGenMETTrg
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
MonitorElement * _meGenMET
void getHLTResults(const edm::TriggerResults &, const edm::TriggerNames &triggerNames)
MonitorElement * _meGenJetPhi
MonitorElement * _meGenHTTrg
MonitorElement * _meRecoJetPhi
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:27
MonitorElement * _meGenJetEtaTrg
MonitorElement * _meGenJetPtTrg
edm::InputTag triggerEventObject_
InputTag of TriggerEventWithRefs to analyze.
TLorentzVector genMet(const HepMC::GenEvent *all, double etamin=-9999., double etamax=9999.)
std::string const & label() const
Definition: InputTag.h:25
MonitorElement * _meRecoHT
MonitorElement * _meGenJetPhiTrgLow
MonitorElement * _meRecoMETTrgLow
MonitorElement * _meGenJetPhiTrg
std::map< std::string, bool >::iterator trig_iter
virtual void analyze(const edm::Event &, const edm::EventSetup &)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:237
MonitorElement * _meRecoJetPt
edm::InputTag CaloJetAlgorithm