CMS 3D CMS Logo

HLTTrackWithHits.h
Go to the documentation of this file.
1 #ifndef HLTrigger_HLTTrackWithHits_H
2 
8 // system include files
9 #include <memory>
10 
11 // user include files
14 
17 
19 
21 
26 
27 class HLTTrackWithHits : public HLTFilter {
28 public:
29  explicit HLTTrackWithHits(const edm::ParameterSet& iConfig) : HLTFilter(iConfig),
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  MinPT_(iConfig.getParameter<double>("MinPT"))
37  {
38  srcToken_ = consumes<reco::TrackCollection>(src_);
39  }
40 
42 
44  {
47  desc.add<edm::InputTag>("src",edm::InputTag(""));
48  desc.add<int>("MinN",0);
49  desc.add<int>("MaxN",99999);
50  desc.add<int>("MinBPX",0);
51  desc.add<int>("MinFPX",0);
52  desc.add<int>("MinPXL",0);
53  desc.add<double>("MinPT",0.);
54  descriptions.add("hltTrackWithHits",desc);
55  }
56 
57 private:
58  virtual bool hltFilter(edm::Event& iEvent, const edm::EventSetup&, trigger::TriggerFilterObjectWithRefs & filterproduct) const override
59  {
61  iEvent.getByToken(srcToken_, oHandle);
62  int s=oHandle->size();
63  int count=0;
64  for (int i=0;i!=s;++i){
65  const reco::Track & track = (*oHandle)[i];
66  if (track.pt() < MinPT_) continue;
67  const reco::HitPattern & hits = track.hitPattern();
68  if ( MinBPX_>0 && hits.numberOfValidPixelBarrelHits() >= MinBPX_ ) { ++count; continue; }
69  if ( MinFPX_>0 && hits.numberOfValidPixelEndcapHits() >= MinFPX_ ) { ++count; continue; }
70  if ( MinPXL_>0 && hits.numberOfValidPixelHits() >= MinPXL_ ) { ++count; continue; }
71  }
72 
73  bool answer=(count>=minN_ && count<=maxN_);
74  LogDebug("HLTTrackWithHits")<<module(iEvent)<<" sees: "<<s<<" objects. Only: "<<count<<" satisfy the hit requirement. Filter answer is: "<<(answer?"true":"false")<<std::endl;
75  return answer;
76  }
77 
81  double MinPT_;
82 };
83 
84 
85 #endif
#define LogDebug(id)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual bool hltFilter(edm::Event &iEvent, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
int iEvent
Definition: GenABIO.cc:230
int numberOfValidPixelBarrelHits() const
Definition: HitPattern.h:843
double pt() const
track transverse momentum
Definition: TrackBase.h:621
int module(edm::Event const &) const
Definition: HLTFilter.cc:52
HLTTrackWithHits(const edm::ParameterSet &iConfig)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:446
void add(std::string const &label, ParameterSetDescription const &psetDescription)
int numberOfValidPixelEndcapHits() const
Definition: HitPattern.h:848
HLT enums.
edm::EDGetTokenT< reco::TrackCollection > srcToken_
int numberOfValidPixelHits() const
Definition: HitPattern.h:838
edm::InputTag src_
answer
Definition: submit.py:44