CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BCToEFilterAlgo.cc
Go to the documentation of this file.
6 
8 
9  //set constants
10  FILTER_ETA_MAX_=2.5;
11  eTThreshold_=(float)iConfig.getParameter<double>("eTThreshold");
12  genParSource_=iC.consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("genParSource"));
13 
14 }
15 
17 }
18 
19 
20 //look for status==1 e coming from b or c hadron
21 //there is an eT threshold on the electron (configurable)
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 }
43 
44 
45 
46 //does this particle have an ancestor (mother, mother of mother, etc.) that is a b or c hadron?
47 //attention: the GenParticle argument must have the motherRef correctly filled for this
48 //to work. That is, you had better have got it out of the genParticles collection
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 }
63 
65  return isBCMeson(gp) || isBCBaryon(gp);
66 }
67 
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 }
79 
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 }
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T getParameter(std::string const &) const
virtual int pdgId() const
PDG identifier.
bool isBCMeson(const reco::GenParticle &gp)
virtual double et() const
transverse energy
bool filter(const edm::Event &iEvent)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
bool isBCBaryon(const reco::GenParticle &gp)
virtual int status() const
status word
bool hasBCAncestors(const reco::GenParticle &gp)
virtual double eta() const
momentum pseudorapidity
BCToEFilterAlgo(const edm::ParameterSet &, edm::ConsumesCollector &&iC)
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::GenParticleCollection > genParSource_
virtual size_t numberOfMothers() const
number of mothers
tuple result
Definition: query.py:137
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
daughters::value_type motherRef(size_type i=0) const
reference to mother at given position
bool isBCHadron(const reco::GenParticle &gp)