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
HFFilter Class Reference

#include <PhysicsTools/HFFilter/src/HFFilter.cc>

Inheritance diagram for HFFilter:
edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

virtual void endJob ()
 
virtual bool filter (edm::Event &, const edm::EventSetup &)
 
 HFFilter (const edm::ParameterSet &)
 
 ~HFFilter ()
 
- Public Member Functions inherited from edm::EDFilter
 EDFilter ()
 
virtual ~EDFilter ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Attributes

double etaMax_
 
edm::InputTag genJetsCollName_
 
double ptMin_
 

Additional Inherited Members

- Public Types inherited from edm::EDFilter
typedef EDFilter ModuleType
 
typedef WorkerT< EDFilterWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDFilter
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDFilter
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

Description: Filter to see if there are heavy flavor GenJets in this event

Implementation: The implementation is simple, it loops over the GenJets and checks if any constituents have a pdg ID that matches a list. It also has a switch to count objects from a gluon parent, so the user can turn off counting gluon splitting.

Definition at line 44 of file HFFilter.h.

Constructor & Destructor Documentation

HFFilter::HFFilter ( const edm::ParameterSet iConfig)
explicit

Definition at line 12 of file HFFilter.cc.

References edm::ParameterSet::getParameter().

13 {
14  genJetsCollName_ = iConfig.getParameter<edm::InputTag>("genJetsCollName");
15  ptMin_ = iConfig.getParameter<double>("ptMin");
16  etaMax_ = iConfig.getParameter<double>("etaMax");
17 }
T getParameter(std::string const &) const
double ptMin_
Definition: HFFilter.h:55
double etaMax_
Definition: HFFilter.h:56
edm::InputTag genJetsCollName_
Definition: HFFilter.h:54
HFFilter::~HFFilter ( )

Definition at line 20 of file HFFilter.cc.

21 {
22 }

Member Function Documentation

void HFFilter::endJob ( void  )
virtual

Reimplemented from edm::EDFilter.

Definition at line 72 of file HFFilter.cc.

72  {
73 }
bool HFFilter::filter ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Implements edm::EDFilter.

Definition at line 32 of file HFFilter.cc.

References JetMCTagUtils::decayFromBHadron(), JetMCTagUtils::decayFromCHadron(), end, edm::Event::getByLabel(), and dt_dqm_sourceclient_common_cff::reco.

33 {
34 
35  // Get the GenJetCollection
36  using namespace edm;
37  using namespace reco;
38  Handle<std::vector<GenJet> > hGenJets;
39  iEvent.getByLabel(genJetsCollName_,hGenJets);
40 
41  // Loop over the GenJetCollection
42  vector<GenJet>::const_iterator ijet = hGenJets->begin();
43  vector<GenJet>::const_iterator end = hGenJets->end();
44  for ( ; ijet != end; ++ijet ) {
45 
46  // Check to make sure the GenJet satisfies kinematic cuts. Ignore those that don't.
47  if ( ijet->pt() < ptMin_ || fabs(ijet->eta()) > etaMax_ ) continue;
48 
49  // Get the constituent particles
50  vector<const GenParticle*> particles = ijet->getGenConstituents ();
51 
52  // Loop over the constituent particles
53  vector<const GenParticle*>::const_iterator genit = particles.begin();
54  vector<const GenParticle*>::const_iterator genend = particles.end();
55  for ( ; genit != genend; ++genit ) {
56 
57  // See if any of them come from B or C hadrons
58  const GenParticle & genitref = **genit;
59  if ( JetMCTagUtils::decayFromBHadron( genitref ) ||
60  JetMCTagUtils::decayFromCHadron( genitref ) ) {
61  return true;
62  }
63  }// end loop over constituents
64  }// end loop over jets
65 
66 
67  return false;
68 }
double ptMin_
Definition: HFFilter.h:55
double etaMax_
Definition: HFFilter.h:56
edm::InputTag genJetsCollName_
Definition: HFFilter.h:54
#define end
Definition: vmac.h:38
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
bool decayFromCHadron(const reco::Candidate &c)
Definition: JetMCTag.cc:61
bool decayFromBHadron(const reco::Candidate &c)
Definition: JetMCTag.cc:40

Member Data Documentation

double HFFilter::etaMax_
private

Definition at line 56 of file HFFilter.h.

edm::InputTag HFFilter::genJetsCollName_
private

Definition at line 54 of file HFFilter.h.

double HFFilter::ptMin_
private

Definition at line 55 of file HFFilter.h.