CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
HeavyQuarkFromMPIFilterAlgo Class Reference

#include <HeavyQuarkFromMPIFilterAlgo.h>

Public Member Functions

bool filter (const edm::Event &iEvent)
 
bool hasMPIAncestor (const reco::GenParticle *)
 
 HeavyQuarkFromMPIFilterAlgo (const edm::ParameterSet &)
 
 ~HeavyQuarkFromMPIFilterAlgo ()
 

Private Attributes

edm::InputTag genParSource_
 
int HeavyQuarkFlavour
 

Detailed Description

Definition at line 23 of file HeavyQuarkFromMPIFilterAlgo.h.

Constructor & Destructor Documentation

HeavyQuarkFromMPIFilterAlgo::HeavyQuarkFromMPIFilterAlgo ( const edm::ParameterSet iConfig)

Definition at line 10 of file HeavyQuarkFromMPIFilterAlgo.cc.

References edm::ParameterSet::getParameter().

10  {
11  HeavyQuarkFlavour=(int)iConfig.getParameter<int>("HQFlavour");
12  genParSource_=iConfig.getParameter<edm::InputTag>("genParSource");
13 }
T getParameter(std::string const &) const
HeavyQuarkFromMPIFilterAlgo::~HeavyQuarkFromMPIFilterAlgo ( )

Definition at line 15 of file HeavyQuarkFromMPIFilterAlgo.cc.

15  {
16 }

Member Function Documentation

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

Definition at line 21 of file HeavyQuarkFromMPIFilterAlgo.cc.

References funct::abs(), edm::RefVector< C, T, F >::begin(), edm::RefVector< C, T, F >::end(), edm::Event::getByLabel(), reco::CompositeRefCandidateT< D >::motherRefVector(), reco::LeafCandidate::pdgId(), and reco::LeafCandidate::status().

21  {
22  bool Veto=false;
24  iEvent.getByLabel(genParSource_,genParsHandle);
25  reco::GenParticleCollection genPars=*genParsHandle;
26 // std::cout<<" in the filter HQF="<<HeavyQuarkFlavour<<std::endl;
27  bool fromMPI=false;
28  for (uint32_t ig=0;ig<genPars.size();ig++) {
29  reco::GenParticle gp=genPars[ig];
30  if (abs(gp.pdgId())==HeavyQuarkFlavour) {
31  if(gp.status()>=30 && gp.status()<40){
32 // cout<<"Found a B with status of 3X (eta="<<gp.eta()<<")==>fromMPI=true"<<endl;
33  fromMPI=true;
34  }
35  if(gp.status()>=40 && fromMPI==false){
36 // cout<<"Found a B with status >= 4X (eta="<<gp.eta()<<")==>Need to check ancestors!"<<endl;
37  const reco::GenParticleRefVector& mothers = gp.motherRefVector();
38 // cout<<"Note: it has "<<mothers.size()<<" mothers"<<endl;
39  for( reco::GenParticleRefVector::const_iterator im = mothers.begin(); im!=mothers.end(); ++im) {
40  const reco::GenParticle& part = **im;
41 // cout<<"--->Going to a mother, having eta="<<part.eta()<<endl;
42  if( hasMPIAncestor( &part) ){
43 // cout<<"------>Found one ancestor with status of 3X (eta="<<gp.eta()<<")==>fromMPI=true"<<endl;
44  fromMPI=true;
45  }
46  }
47  }
48  if(fromMPI)Veto=true;
49  else Veto=false;
50  }
51  }
52 // cout<<"RETURN "<<Veto<<endl;
53  return Veto;
54 }
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
virtual int pdgId() const
PDG identifier.
virtual int status() const
status word
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:255
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:250
bool hasMPIAncestor(const reco::GenParticle *)
const mothers & motherRefVector() const
references to mothers
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
part
Definition: HCALResponse.h:20
bool HeavyQuarkFromMPIFilterAlgo::hasMPIAncestor ( const reco::GenParticle particle)

Definition at line 56 of file HeavyQuarkFromMPIFilterAlgo.cc.

References edm::RefVector< C, T, F >::begin(), edm::RefVector< C, T, F >::end(), reco::CompositeRefCandidateT< D >::motherRefVector(), and reco::LeafCandidate::status().

56  {
57  if( particle->status() >=30 && particle->status()<40 ){
58 // cout<<"------->Mother found with eta="<<particle->eta()<<" and status "<<particle->status()<<"==> returning true!"<<endl;
59  return true;
60  }
61 // cout<<"------->in hasMPIAncestor, current particle has eta="<<particle->eta()<<" and status!=3X ==> checking next mothers"<<endl;
62  const reco::GenParticleRefVector& mothers = particle->motherRefVector();
63 // cout<<"------>Note: it has "<<mothers.size()<<" mothers"<<endl;
64  for( reco::GenParticleRefVector::const_iterator im = mothers.begin(); im!=mothers.end(); ++im) {
65  const reco::GenParticle& part = **im;
66  if( hasMPIAncestor( &part )){
67 // cout<<"--------->Found one ancestor with status of 3X ==>returning true"<<endl;
68  return true;
69  }
70  }
71 // cout<<"------>No 30's ancestor found ==>returning false"<<endl;
72  return false;
73 }
virtual int status() const
status word
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:255
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:250
bool hasMPIAncestor(const reco::GenParticle *)
const mothers & motherRefVector() const
references to mothers
part
Definition: HCALResponse.h:20

Member Data Documentation

edm::InputTag HeavyQuarkFromMPIFilterAlgo::genParSource_
private

Definition at line 33 of file HeavyQuarkFromMPIFilterAlgo.h.

int HeavyQuarkFromMPIFilterAlgo::HeavyQuarkFlavour
private

Definition at line 32 of file HeavyQuarkFromMPIFilterAlgo.h.