CMS 3D CMS Logo

TopGenEventAnalyzer.cc
Go to the documentation of this file.
2 
4  inputGenEventToken_(consumes<TtGenEvent>(cfg.getParameter<edm::InputTag>("genEvent")))
5 {
7  nLep_ = fs->make<TH1F>("nLep", "N(Lepton)", 5, 0., 5.);
8  topPt_ = fs->make<TH1F>("topPt", "pt (top)", 100, 0., 500.);
9  topEta_ = fs->make<TH1F>("topEta", "eta(top)", 40, -5., 5.);
10  topPhi_ = fs->make<TH1F>("topPhi", "phi(top)", 60, -3.5, 3.5);
11  topBarPt_ = fs->make<TH1F>("topBarPt", "pt (topBar)", 100, 0., 500.);
12  topBarEta_ = fs->make<TH1F>("topBarEta", "eta(topBar)", 40, -5., 5.);
13  topBarPhi_ = fs->make<TH1F>("topBarPhi", "phi(topBar)", 60, -3.5, 3.5);
14  ttbarPt_ = fs->make<TH1F>("ttbarPt", "pt (ttbar)", 100, 0., 500.);
15  ttbarEta_ = fs->make<TH1F>("ttbarEta", "eta(ttbar)", 40, -5., 5.);
16  ttbarPhi_ = fs->make<TH1F>("ttbarPhi", "phi(ttbar)", 60, -3.5, 3.5);
17  prodChan_ = fs->make<TH1F>("prodChan", "production mode", 3, 0, 3);
18  prodChan_->GetXaxis()->SetBinLabel(1, "gg" );
19  prodChan_->GetXaxis()->SetBinLabel(2, "qqbar");
20  prodChan_->GetXaxis()->SetBinLabel(3, "other");
21 }
22 
24 {
25 }
26 
27 void
29 {
31  evt.getByToken(inputGenEventToken_, genEvent);
32 
33  if(!genEvent->isTtBar())
34  return;
35 
36  if(genEvent->fromGluonFusion())
37  prodChan_->Fill("gg", 1);
38  else if(genEvent->fromQuarkAnnihilation())
39  prodChan_->Fill("qqbar", 1);
40  else
41  prodChan_->Fill("other", 1);
42 
43  // fill BR's
44  nLep_ ->Fill(genEvent->numberOfLeptons());
45 
46  //fill top kinematic
47  topPt_ ->Fill(genEvent->top ()->pt ());
48  topEta_ ->Fill(genEvent->top ()->eta());
49  topPhi_ ->Fill(genEvent->top ()->phi());
50  topBarPt_ ->Fill(genEvent->topBar()->pt ());
51  topBarEta_->Fill(genEvent->topBar()->eta());
52  topBarPhi_->Fill(genEvent->topBar()->phi());
53 
54  //fill ttbar kinematics
55  ttbarPt_ ->Fill(genEvent->topPair()->pt() );
56  ttbarEta_->Fill(genEvent->topPair()->eta());
57  ttbarPhi_->Fill(genEvent->topPair()->phi());
58 }
59 
61 {
62 }
63 
65 {
66 }
double eta() const final
momentum pseudorapidity
TopGenEventAnalyzer(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
bool fromQuarkAnnihilation() const
check if the tops were produced from qqbar
Definition: TtGenEvent.cc:29
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
double pt() const final
transverse momentum
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
const reco::GenParticle * top() const
return top if available; 0 else
Definition: TopGenEvent.h:104
Class derived from the TopGenEvent for ttbar events.
Definition: TtGenEvent.h:18
bool fromGluonFusion() const
check if the tops were produced from a pair of gluons
Definition: TtGenEvent.cc:19
void analyze(const edm::Event &, const edm::EventSetup &) override
void beginJob() override
edm::EDGetTokenT< TtGenEvent > inputGenEventToken_
HLT enums.
const math::XYZTLorentzVector * topPair() const
return combined 4-vector of top and topBar
Definition: TtGenEvent.h:87
double phi() const final
momentum azimuthal angle
const reco::GenParticle * topBar() const
return anti-top if available; 0 else
Definition: TopGenEvent.h:106
bool isTtBar() const
check if the event can be classified as ttbar
Definition: TtGenEvent.h:30
int numberOfLeptons(bool fromWBoson=true) const
return number of leptons in the decay chain
Definition: TopGenEvent.cc:48