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