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 Member Functions | Private Attributes
HLTTrackWithHits Class Reference

#include <HLTTrackWithHits.h>

Inheritance diagram for HLTTrackWithHits:
HLTFilter edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

 HLTTrackWithHits (const edm::ParameterSet &iConfig)
 
 ~HLTTrackWithHits ()
 
- Public Member Functions inherited from HLTFilter
 HLTFilter ()
 
virtual ~HLTFilter ()
 
- 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 Member Functions

virtual void endJob ()
 
virtual bool filter (edm::Event &iEvent, const edm::EventSetup &)
 

Private Attributes

int maxN_
 
int MinBPX_
 
int MinFPX_
 
int minN_
 
int MinPXL_
 
edm::InputTag src_
 

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)
 
- 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: templated EDFilter to count the number of tracks with a given hit requirement

Author
Jean-Roch Vlimant

Definition at line 27 of file HLTTrackWithHits.h.

Constructor & Destructor Documentation

HLTTrackWithHits::HLTTrackWithHits ( const edm::ParameterSet iConfig)
inlineexplicit

Definition at line 29 of file HLTTrackWithHits.h.

29  :
30  src_(iConfig.getParameter<edm::InputTag>("src")),
31  minN_(iConfig.getParameter<int>("MinN")),
32  maxN_(iConfig.getParameter<int>("MaxN")),
33  MinBPX_(iConfig.getParameter<int>("MinBPX")),
34  MinFPX_(iConfig.getParameter<int>("MinFPX")),
35  MinPXL_(iConfig.getParameter<int>("MinPXL"))
36  {
37  produces<trigger::TriggerFilterObjectWithRefs>();
38  };
T getParameter(std::string const &) const
edm::InputTag src_
HLTTrackWithHits::~HLTTrackWithHits ( )
inline

Definition at line 40 of file HLTTrackWithHits.h.

40 {};

Member Function Documentation

virtual void HLTTrackWithHits::endJob ( void  )
inlineprivatevirtual

Reimplemented from edm::EDFilter.

Definition at line 66 of file HLTTrackWithHits.h.

66 {};
virtual bool HLTTrackWithHits::filter ( edm::Event iEvent,
const edm::EventSetup  
)
inlineprivatevirtual

Implements HLTFilter.

Definition at line 43 of file HLTTrackWithHits.h.

References submit::answer, prof2calltree::count, edm::Event::getByLabel(), reco::TrackBase::hitPattern(), i, LogDebug, maxN_, MinBPX_, MinFPX_, minN_, MinPXL_, module(), reco::HitPattern::numberOfValidPixelBarrelHits(), reco::HitPattern::numberOfValidPixelEndcapHits(), reco::HitPattern::numberOfValidPixelHits(), path(), edm::Event::put(), asciidump::s, and src_.

44  {
45  // The filtered object. which is put empty.
46  std::auto_ptr<trigger::TriggerFilterObjectWithRefs> filterproduct (new trigger::TriggerFilterObjectWithRefs(path(),module()));
47 
49  iEvent.getByLabel(src_, oHandle);
50  int s=oHandle->size();
51  int count=0;
52  for (int i=0;i!=s;++i){
53  const reco::Track & track = (*oHandle)[i];
54  const reco::HitPattern & hits = track.hitPattern();
55  if ( MinBPX_>0 && hits.numberOfValidPixelBarrelHits() >= MinBPX_ ) count++; continue;
56  if ( MinFPX_>0 && hits.numberOfValidPixelEndcapHits() >= MinFPX_ ) count++; continue;
57  if ( MinPXL_>0 && hits.numberOfValidPixelHits() >= MinPXL_ ) count++; continue;
58  }
59 
60  bool answer=(count>=minN_ && count<=maxN_);
61  LogDebug("HLTTrackWithHits")<<module()<<" sees: "<<s<<" objects. Only: "<<count<<" satisfy the hit requirement. Filter answer is: "<<(answer?"true":"false")<<std::endl;
62 
63  iEvent.put(filterproduct);
64  return answer;
65  }
#define LogDebug(id)
answer
Definition: submit.py:44
int i
Definition: DBlmapReader.cc:9
int module() const
Definition: HLTadd.h:12
int path() const
Definition: HLTadd.h:3
int numberOfValidPixelEndcapHits() const
Definition: HitPattern.cc:385
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
int numberOfValidPixelHits() const
Definition: HitPattern.cc:359
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:223
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
int numberOfValidPixelBarrelHits() const
Definition: HitPattern.cc:372
string s
Definition: asciidump.py:422
edm::InputTag src_

Member Data Documentation

int HLTTrackWithHits::maxN_
private

Definition at line 69 of file HLTTrackWithHits.h.

Referenced by filter().

int HLTTrackWithHits::MinBPX_
private

Definition at line 69 of file HLTTrackWithHits.h.

Referenced by filter().

int HLTTrackWithHits::MinFPX_
private

Definition at line 69 of file HLTTrackWithHits.h.

Referenced by filter().

int HLTTrackWithHits::minN_
private

Definition at line 69 of file HLTTrackWithHits.h.

Referenced by filter().

int HLTTrackWithHits::MinPXL_
private

Definition at line 69 of file HLTTrackWithHits.h.

Referenced by filter().

edm::InputTag HLTTrackWithHits::src_
private

Definition at line 66 of file HLTTrackWithHits.h.

Referenced by filter().