CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
TopGenEventAnalyzer Class Reference

#include <TopGenEventAnalyzer.h>

Inheritance diagram for TopGenEventAnalyzer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 TopGenEventAnalyzer (const edm::ParameterSet &)
 
 ~TopGenEventAnalyzer () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

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

Private Attributes

edm::EDGetTokenT< TtGenEventinputGenEventToken_
 
TH1F * nLep_
 
TH1F * prodChan_
 
TH1F * topBarEta_
 
TH1F * topBarPhi_
 
TH1F * topBarPt_
 
TH1F * topEta_
 
TH1F * topPhi_
 
TH1F * topPt_
 
TH1F * ttbarEta_
 
TH1F * ttbarPhi_
 
TH1F * ttbarPt_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 16 of file TopGenEventAnalyzer.h.

Constructor & Destructor Documentation

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

Definition at line 3 of file TopGenEventAnalyzer.cc.

References TFileService::make(), nLep_, prodChan_, topBarEta_, topBarPhi_, topBarPt_, topEta_, topPhi_, topPt_, ttbarEta_, ttbarPhi_, and ttbarPt_.

3  :
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 }
T getParameter(std::string const &) const
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
edm::EDGetTokenT< TtGenEvent > inputGenEventToken_
TopGenEventAnalyzer::~TopGenEventAnalyzer ( )
override

Definition at line 23 of file TopGenEventAnalyzer.cc.

24 {
25 }

Member Function Documentation

void TopGenEventAnalyzer::analyze ( const edm::Event evt,
const edm::EventSetup setup 
)
overrideprivate

Definition at line 28 of file TopGenEventAnalyzer.cc.

References reco::LeafCandidate::eta(), TtGenEvent::fromGluonFusion(), TtGenEvent::fromQuarkAnnihilation(), nano_cff::genEvent, edm::Event::getByToken(), inputGenEventToken_, TtGenEvent::isTtBar(), nLep_, TopGenEvent::numberOfLeptons(), reco::LeafCandidate::phi(), prodChan_, reco::LeafCandidate::pt(), TopGenEvent::top(), TopGenEvent::topBar(), topBarEta_, topBarPhi_, topBarPt_, topEta_, TtGenEvent::topPair(), topPhi_, topPt_, ttbarEta_, ttbarPhi_, and ttbarPt_.

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 }
double eta() const final
momentum pseudorapidity
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
bool fromQuarkAnnihilation() const
check if the tops were produced from qqbar
Definition: TtGenEvent.cc:29
double pt() const final
transverse momentum
const reco::GenParticle * top() const
return top if available; 0 else
Definition: TopGenEvent.h:104
bool fromGluonFusion() const
check if the tops were produced from a pair of gluons
Definition: TtGenEvent.cc:19
edm::EDGetTokenT< TtGenEvent > inputGenEventToken_
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
void TopGenEventAnalyzer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 60 of file TopGenEventAnalyzer.cc.

61 {
62 }
void TopGenEventAnalyzer::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 64 of file TopGenEventAnalyzer.cc.

65 {
66 }

Member Data Documentation

edm::EDGetTokenT<TtGenEvent> TopGenEventAnalyzer::inputGenEventToken_
private

Definition at line 29 of file TopGenEventAnalyzer.h.

Referenced by analyze().

TH1F* TopGenEventAnalyzer::nLep_
private

Definition at line 31 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::prodChan_
private

Definition at line 41 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topBarEta_
private

Definition at line 36 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topBarPhi_
private

Definition at line 37 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topBarPt_
private

Definition at line 35 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topEta_
private

Definition at line 33 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topPhi_
private

Definition at line 34 of file TopGenEventAnalyzer.h.

Referenced by analyze(), and TopGenEventAnalyzer().

TH1F* TopGenEventAnalyzer::topPt_
private

Definition at line 32 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::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().