CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

TopGenEventAnalyzer Class Reference

#include <TopGenEventAnalyzer.h>

Inheritance diagram for TopGenEventAnalyzer:
edm::EDAnalyzer

List of all members.

Public Member Functions

 TopGenEventAnalyzer (const edm::ParameterSet &)
 ~TopGenEventAnalyzer ()

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob ()
virtual void endJob ()

Private Attributes

edm::InputTag inputGenEvent_
TH1F * nLep_
TH1F * topBarEta_
TH1F * topBarPhi_
TH1F * topBarPt_
TH1F * topEta_
TH1F * topPhi_
TH1F * topPt_
TH1F * ttbarEta_
TH1F * ttbarPhi_
TH1F * ttbarPt_

Detailed Description

Definition at line 14 of file TopGenEventAnalyzer.h.


Constructor & Destructor Documentation

TopGenEventAnalyzer::TopGenEventAnalyzer ( const edm::ParameterSet cfg) [explicit]

Definition at line 4 of file TopGenEventAnalyzer.cc.

References nLep_, topBarEta_, topBarPhi_, topBarPt_, topEta_, topPhi_, topPt_, ttbarEta_, ttbarPhi_, and ttbarPt_.

                                                                  :
  inputGenEvent_(cfg.getParameter<edm::InputTag>("genEvent"))
{ 
  edm::Service<TFileService> fs;
  nLep_      = fs->make<TH1F>("nLep",      "N(Lepton)",     5,   0.,   5.);
  topPt_     = fs->make<TH1F>("topPt",     "pt (top)",    100,   0., 500.);
  topEta_    = fs->make<TH1F>("topEta",    "eta(top)",     40,  -5.,   5.);
  topPhi_    = fs->make<TH1F>("topPhi",    "phi(top)",     60, -3.5,  3.5);
  topBarPt_  = fs->make<TH1F>("topBarPt",  "pt (topBar)", 100,   0., 500.);
  topBarEta_ = fs->make<TH1F>("topBarEta", "eta(topBar)",  40,  -5.,   5.);
  topBarPhi_ = fs->make<TH1F>("topBarPhi", "phi(topBar)",  60, -3.5,  3.5);
  ttbarPt_   = fs->make<TH1F>("ttbarPt",   "pt (ttbar)",  100,   0., 500.);
  ttbarEta_  = fs->make<TH1F>("ttbarEta",  "eta(ttbar)",   40,  -5.,   5.);
  ttbarPhi_  = fs->make<TH1F>("ttbarPhi",  "phi(ttbar)",   60, -3.5,  3.5);
}
TopGenEventAnalyzer::~TopGenEventAnalyzer ( )

Definition at line 20 of file TopGenEventAnalyzer.cc.

{
}

Member Function Documentation

void TopGenEventAnalyzer::analyze ( const edm::Event evt,
const edm::EventSetup setup 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 25 of file TopGenEventAnalyzer.cc.

References MCTruth::genEvent, edm::Event::getByLabel(), inputGenEvent_, nLep_, p4, topBarEta_, topBarPhi_, topBarPt_, topEta_, topPhi_, topPt_, ttbarEta_, ttbarPhi_, and ttbarPt_.

{
  edm::Handle<TtGenEvent> genEvent;
  evt.getByLabel(inputGenEvent_, genEvent);

  // fill BR's
  nLep_  ->Fill(genEvent->numberOfLeptons());

  //fill top kinematic
  topPt_    ->Fill(genEvent->top   ()->pt ());
  topEta_   ->Fill(genEvent->top   ()->eta());
  topPhi_   ->Fill(genEvent->top   ()->phi());
  topBarPt_ ->Fill(genEvent->topBar()->pt ());
  topBarEta_->Fill(genEvent->topBar()->eta());
  topBarPhi_->Fill(genEvent->topBar()->phi());

  //fill ttbar kinematics
  reco::Particle::LorentzVector p4 = genEvent->top()->p4()+genEvent->topBar()->p4();
  ttbarPt_ ->Fill(p4.pt() );
  ttbarEta_->Fill(p4.eta());
  ttbarPhi_->Fill(p4.phi());
}
void TopGenEventAnalyzer::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 48 of file TopGenEventAnalyzer.cc.

{  
} 
void TopGenEventAnalyzer::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 52 of file TopGenEventAnalyzer.cc.

{
}

Member Data Documentation

Definition at line 27 of file TopGenEventAnalyzer.h.

Referenced by analyze().

TH1F* TopGenEventAnalyzer::nLep_ [private]

Definition at line 29 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

Definition at line 34 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

Definition at line 35 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

Definition at line 33 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

Definition at line 31 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

Definition at line 32 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topPt_ [private]

Definition at line 30 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

Definition at line 37 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

Definition at line 38 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

Definition at line 36 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().