CMS 3D CMS Logo

TopGenEventAnalyzer Class Reference

#include <TopQuarkAnalysis/Examples/plugins/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 (const edm::EventSetup &)
virtual void endJob ()

Private Attributes

edm::InputTag inputGenEvent_
TH1F * nLep_
TH1F * topBarEta_
TH1F * topBarMass_
TH1F * topBarPhi_
TH1F * topBarPt_
TH1F * topEta_
TH1F * topMass_
TH1F * topPhi_
TH1F * topPt_
TH1F * ttbarEta_
TH1F * ttbarMass_
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 6 of file TopGenEventAnalyzer.cc.

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

00006                                                                   :
00007   inputGenEvent_(cfg.getParameter<edm::InputTag>("genEvent"))
00008 { 
00009   edm::Service<TFileService> fs;
00010   nLep_      = fs->make<TH1F>("nLep",      "N(Lepton)",     5,   0.,   5.);
00011   topPt_     = fs->make<TH1F>("topPt",     "pt (top)",    100,   0., 500.);
00012   topEta_    = fs->make<TH1F>("topEta",    "eta(top)",     40,  -5.,   5.);
00013   topPhi_    = fs->make<TH1F>("topPhi",    "phi(top)",     60, -3.5,  3.5);
00014   topMass_   = fs->make<TH1F>("topMass",   "mass(top)",   150, 100., 250.);
00015   topBarPt_  = fs->make<TH1F>("topBarPt",  "pt (topBar)", 100,   0., 500.);
00016   topBarEta_ = fs->make<TH1F>("topBarEta", "eta(topBar)",  40,  -5.,   5.);
00017   topBarPhi_ = fs->make<TH1F>("topBarPhi", "phi(topBar)",  60, -3.5,  3.5);
00018   topBarMass_= fs->make<TH1F>("topBarMass","mass(top)",   150, 100., 250.);
00019   ttbarPt_   = fs->make<TH1F>("ttbarPt",   "pt (ttbar)",  100,   0., 500.);
00020   ttbarEta_  = fs->make<TH1F>("ttbarEta",  "eta(ttbar)",   40,  -5.,   5.);
00021   ttbarPhi_  = fs->make<TH1F>("ttbarPhi",  "phi(ttbar)",   60, -3.5,  3.5);
00022   ttbarMass_ = fs->make<TH1F>("ttbarMass", "mass(ttbar)", 150, 100., 250.);
00023 }

TopGenEventAnalyzer::~TopGenEventAnalyzer (  ) 

Definition at line 25 of file TopGenEventAnalyzer.cc.

00026 {
00027 }


Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 30 of file TopGenEventAnalyzer.cc.

References MCTruth2::genEvent, edm::Event::getByLabel(), inputGenEvent_, funct::log(), nLep_, p4, topBarEta_, topBarMass_, topBarPhi_, topBarPt_, topEta_, topMass_, topPhi_, topPt_, ttbarEta_, ttbarMass_, ttbarPhi_, and ttbarPt_.

00031 {
00032   edm::Handle<TtGenEvent> genEvent;
00033   evt.getByLabel(inputGenEvent_, genEvent);
00034 
00035   if(genEvent->isFullLeptonic(true)){
00036     genEvent->dumpEventContent();
00037     edm::LogVerbatim log("TopGenEventAnalyzer::selection");
00038     log << "!!! - is full-leptonic - !!! \n";
00039   }
00040 
00041   // fill BR's
00042   nLep_  ->Fill(genEvent->numberOfLeptons());
00043 
00044   //fill top kinematic
00045   topPt_     ->Fill(genEvent->top   ()->pt  ());
00046   topEta_    ->Fill(genEvent->top   ()->eta ());
00047   topPhi_    ->Fill(genEvent->top   ()->phi ());
00048   topMass_   ->Fill(genEvent->top   ()->mass());
00049   topBarPt_  ->Fill(genEvent->topBar()->pt  ());
00050   topBarEta_ ->Fill(genEvent->topBar()->eta ());
00051   topBarPhi_ ->Fill(genEvent->topBar()->phi ());
00052   topBarMass_->Fill(genEvent->topBar()->mass());
00053 
00054   //fill ttbar kinematics
00055   reco::Particle::LorentzVector p4 = genEvent->top()->p4()+genEvent->topBar()->p4();
00056   ttbarPt_  ->Fill(p4.pt  ());
00057   ttbarEta_ ->Fill(p4.eta ());
00058   ttbarPhi_ ->Fill(p4.phi ());
00059   ttbarMass_->Fill(p4.mass());
00060 }

void TopGenEventAnalyzer::beginJob ( const edm::EventSetup  )  [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 62 of file TopGenEventAnalyzer.cc.

00063 {  
00064 } 

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

Reimplemented from edm::EDAnalyzer.

Definition at line 66 of file TopGenEventAnalyzer.cc.

00067 {
00068 }


Member Data Documentation

edm::InputTag TopGenEventAnalyzer::inputGenEvent_ [private]

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().

TH1F* TopGenEventAnalyzer::topBarEta_ [private]

Definition at line 35 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topBarMass_ [private]

Definition at line 37 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topBarPhi_ [private]

Definition at line 36 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topBarPt_ [private]

Definition at line 34 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topEta_ [private]

Definition at line 31 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topMass_ [private]

Definition at line 33 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topPhi_ [private]

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().

TH1F* TopGenEventAnalyzer::ttbarEta_ [private]

Definition at line 39 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::ttbarMass_ [private]

Definition at line 41 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::ttbarPhi_ [private]

Definition at line 40 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::ttbarPt_ [private]

Definition at line 38 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:33:41 2009 for CMSSW by  doxygen 1.5.4