CMS 3D CMS Logo

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

#include <BCToEFilterAlgo.h>

Public Member Functions

 BCToEFilterAlgo (const edm::ParameterSet &, edm::ConsumesCollector &&iC)
 
bool filter (const edm::Event &iEvent)
 
bool hasBCAncestors (const reco::GenParticle &gp)
 
 ~BCToEFilterAlgo ()
 

Private Member Functions

bool isBCBaryon (const reco::GenParticle &gp)
 
bool isBCHadron (const reco::GenParticle &gp)
 
bool isBCMeson (const reco::GenParticle &gp)
 

Private Attributes

float eTThreshold_
 
float FILTER_ETA_MAX_
 
edm::EDGetTokenT< reco::GenParticleCollectiongenParSource_
 

Detailed Description

BCToEFilterAlgo returns true for events that have an electron, above configurable eT threshold and within |eta|<2.5, that has an ancestor of a b or c quark

Author
J Lamb, UCSB

Definition at line 26 of file BCToEFilterAlgo.h.

Constructor & Destructor Documentation

BCToEFilterAlgo::BCToEFilterAlgo ( const edm::ParameterSet iConfig,
edm::ConsumesCollector &&  iC 
)

Definition at line 7 of file BCToEFilterAlgo.cc.

References eTThreshold_, FILTER_ETA_MAX_, objects.autophobj::float, genParSource_, and edm::ParameterSet::getParameter().

7  {
8 
9  //set constants
10  FILTER_ETA_MAX_=2.5;
11  eTThreshold_=(float)iConfig.getParameter<double>("eTThreshold");
13 
14 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::GenParticleCollection > genParSource_
BCToEFilterAlgo::~BCToEFilterAlgo ( )

Definition at line 16 of file BCToEFilterAlgo.cc.

16  {
17 }

Member Function Documentation

bool BCToEFilterAlgo::filter ( const edm::Event iEvent)

Definition at line 22 of file BCToEFilterAlgo.cc.

References funct::abs(), reco::LeafCandidate::et(), reco::LeafCandidate::eta(), eTThreshold_, FILTER_ETA_MAX_, genParSource_, edm::Event::getByToken(), runTauDisplay::gp, hasBCAncestors(), reco::LeafCandidate::pdgId(), mps_fire::result, and reco::LeafCandidate::status().

22  {
23 
24  bool result=false;
25 
26 
27 
29  iEvent.getByToken(genParSource_,genParsHandle);
30  reco::GenParticleCollection genPars=*genParsHandle;
31 
32  for (uint32_t ig=0;ig<genPars.size();ig++) {
33  reco::GenParticle gp=genPars.at(ig);
34  if (gp.status()==1 && abs(gp.pdgId())==11 && gp.et()>eTThreshold_ && fabs(gp.eta())<FILTER_ETA_MAX_) {
35  if (hasBCAncestors(gp)) {
36  result=true;
37  }
38  }
39  }
40  return result;
41 
42 }
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
int pdgId() const final
PDG identifier.
double eta() const final
momentum pseudorapidity
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
bool hasBCAncestors(const reco::GenParticle &gp)
edm::EDGetTokenT< reco::GenParticleCollection > genParSource_
double et() const final
transverse energy
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int status() const final
status word
bool BCToEFilterAlgo::hasBCAncestors ( const reco::GenParticle gp)

Definition at line 49 of file BCToEFilterAlgo.cc.

References isBCHadron(), reco::CompositeRefCandidateT< D >::motherRef(), and reco::CompositeRefCandidateT< D >::numberOfMothers().

Referenced by filter().

49  {
50 
51  //stopping condition: this particle is a b or c hadron
52  if (isBCHadron(gp)) return true;
53  //stopping condition: this particle has no mothers
54  if (gp.numberOfMothers()==0) return false;
55  //otherwise continue
56  bool retval=false;
57  for (uint32_t im=0;im<gp.numberOfMothers();im++) {
58  retval=retval || hasBCAncestors(*gp.motherRef(im));
59  }
60  return retval;
61 
62 }
size_t numberOfMothers() const override
number of mothers
bool hasBCAncestors(const reco::GenParticle &gp)
daughters::value_type motherRef(size_type i=0) const
reference to mother at given position
bool isBCHadron(const reco::GenParticle &gp)
bool BCToEFilterAlgo::isBCBaryon ( const reco::GenParticle gp)
private

Definition at line 80 of file BCToEFilterAlgo.cc.

References funct::abs(), BPhysicsValidation_cfi::pdgid, and reco::LeafCandidate::pdgId().

Referenced by isBCHadron().

80  {
81 
82  uint32_t pdgid=abs(gp.pdgId());
83  uint32_t thousands=pdgid%10000;
84  if (thousands>=4000 && thousands <6000) {
85  return true;
86  } else {
87  return false;
88  }
89 
90 }
int pdgId() const final
PDG identifier.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool BCToEFilterAlgo::isBCHadron ( const reco::GenParticle gp)
private

Definition at line 64 of file BCToEFilterAlgo.cc.

References isBCBaryon(), and isBCMeson().

Referenced by hasBCAncestors().

64  {
65  return isBCMeson(gp) || isBCBaryon(gp);
66 }
bool isBCMeson(const reco::GenParticle &gp)
bool isBCBaryon(const reco::GenParticle &gp)
bool BCToEFilterAlgo::isBCMeson ( const reco::GenParticle gp)
private

Definition at line 68 of file BCToEFilterAlgo.cc.

References funct::abs(), BPhysicsValidation_cfi::pdgid, and reco::LeafCandidate::pdgId().

Referenced by isBCHadron().

68  {
69 
70  uint32_t pdgid=abs(gp.pdgId());
71  uint32_t hundreds=pdgid%1000;
72  if (hundreds>=400 && hundreds<600) {
73  return true;
74  } else {
75  return false;
76  }
77 
78 }
int pdgId() const final
PDG identifier.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

Member Data Documentation

float BCToEFilterAlgo::eTThreshold_
private

Definition at line 47 of file BCToEFilterAlgo.h.

Referenced by BCToEFilterAlgo(), and filter().

float BCToEFilterAlgo::FILTER_ETA_MAX_
private

Definition at line 45 of file BCToEFilterAlgo.h.

Referenced by BCToEFilterAlgo(), and filter().

edm::EDGetTokenT<reco::GenParticleCollection> BCToEFilterAlgo::genParSource_
private

Definition at line 48 of file BCToEFilterAlgo.h.

Referenced by BCToEFilterAlgo(), and filter().