CMS 3D CMS Logo

TTbar_GenJetAnalyzer.cc
Go to the documentation of this file.
3 
5  jets_(iConfig.getParameter<edm::InputTag>("jets")),
6  genEventInfoProductTag_(iConfig.getParameter<edm::InputTag>("genEventInfoProductTag"))
7 {
8 
9  genEventInfoProductTagToken_=consumes<GenEventInfoProduct>(genEventInfoProductTag_);
10  jetsToken_=consumes<std::vector<reco::GenJet> >(jets_);
11 
12 }
13 
14 
16 
17 
19 {
20  using namespace edm;
21 
22  // --- the MC weights ---
24  iEvent.getByToken(genEventInfoProductTagToken_, evt_info);
25  if(!evt_info.isValid()) return;
26  weight = evt_info->weight() ;
27 
28  // Gather information in the GenJet collection
30  iEvent.getByToken(jetsToken_, jets );
31 
32  if(!jets.isValid()) return;
33  // loop Jet collection and fill histograms
34  int njets = 0;
35  for(std::vector<reco::GenJet>::const_iterator jet_it=jets->begin(); jet_it!=jets->end(); ++jet_it){
36 
37  ++njets;
38 
39  hists_["jetPtAll" ]->Fill( jet_it->pt() , weight );
40  hists_["jetEtaAll"]->Fill( jet_it->eta(), weight );
41 
42  if (njets == 1) { hists_["jetPt1" ]->Fill( jet_it->pt() , weight );
43  hists_["jetEta1"]->Fill( jet_it->eta(), weight );
44  }
45  if (njets == 2) { hists_["jetPt2" ]->Fill( jet_it->pt() , weight );
46  hists_["jetEta2"]->Fill( jet_it->eta(), weight );
47  }
48  if (njets == 3) { hists_["jetPt3" ]->Fill( jet_it->pt() , weight );
49  hists_["jetEta3"]->Fill( jet_it->eta(), weight );
50  }
51  if (njets == 4) { hists_["jetPt4" ]->Fill( jet_it->pt() , weight );
52  hists_["jetEta4"]->Fill( jet_it->eta(), weight );
53  }
54  }
55 
56 }
57 
58 
60  DQMHelper dqm(&i); i.setCurrentFolder("Generator/TTbar");
61  hists_["jetPtAll" ] = dqm.book1dHisto("TTbar_jetPtAll" , "pt" , 1000, 0., 1000.,"P_{t}^{All-Jets} (GeV)","Number of Events");
62  hists_["jetPt1" ] = dqm.book1dHisto("TTbar_jetPt1" , "pt" , 1000, 0., 1000.,"P_{t}^{1st-Jet} (GeV)","Number of Events");
63  hists_["jetPt2" ] = dqm.book1dHisto("TTbar_jetPt2" , "pt" , 1000, 0., 1000.,"P_{t}^{2nd-Jet} (GeV)","Number of Events");
64  hists_["jetPt3" ] = dqm.book1dHisto("TTbar_jetPt3" , "pt" , 1000, 0., 1000.,"P_{t}^{3rd-Jet} (GeV)","Number of Events");
65  hists_["jetPt4" ] = dqm.book1dHisto("TTbar_jetPt4" , "pt" , 1000, 0., 1000.,"P_{t}^{4th-Jet} (GeV)","Number of Events");
66 
67  hists_["jetEtaAll"] = dqm.book1dHisto("TTbar_jetEtaAll", "eta", 100, -5., 5.,"#eta^{All-Jets}","Number of Events");
68  hists_["jetEta1" ] = dqm.book1dHisto("TTbar_jetEta1" , "eta", 100, -5., 5.,"#eta^{1st-Jet}","Number of Events");
69  hists_["jetEta2" ] = dqm.book1dHisto("TTbar_jetEta2" , "eta", 100, -5., 5.,"#eta^{2nd-Jet}","Number of Events");
70  hists_["jetEta3" ] = dqm.book1dHisto("TTbar_jetEta3" , "eta", 100, -5., 5.,"#eta^{3rd-Jet}","Number of Events");
71  hists_["jetEta4" ] = dqm.book1dHisto("TTbar_jetEta4" , "eta", 100, -5., 5.,"#eta^{4th-Jet}","Number of Events");
72 }
73 
edm::EDGetTokenT< std::vector< reco::GenJet > > jetsToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
virtual void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
Definition: weight.py:1
double weight() const
MonitorElement * book1dHisto(std::string name, std::string title, int n, double xmin, double xmax, std::string xaxis, std::string yaxis)
Definition: DQMHelper.cc:12
int iEvent
Definition: GenABIO.cc:230
vector< PseudoJet > jets
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
edm::InputTag genEventInfoProductTag_
edm::EDGetTokenT< GenEventInfoProduct > genEventInfoProductTagToken_
bool isValid() const
Definition: HandleBase.h:74
std::map< std::string, MonitorElement * > hists_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
HLT enums.
Definition: Run.h:42
TTbar_GenJetAnalyzer(const edm::ParameterSet &)