CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
24 
25 class HLTTrackWithHits : public HLTFilter {
26 public:
28  : HLTFilter(iConfig),
29  src_(iConfig.getParameter<edm::InputTag>("src")),
30  minN_(iConfig.getParameter<int>("MinN")),
31  maxN_(iConfig.getParameter<int>("MaxN")),
32  MinBPX_(iConfig.getParameter<int>("MinBPX")),
33  MinFPX_(iConfig.getParameter<int>("MinFPX")),
34  MinPXL_(iConfig.getParameter<int>("MinPXL")),
35  MinPT_(iConfig.getParameter<double>("MinPT")) {
36  srcToken_ = consumes<reco::TrackCollection>(src_);
37  }
38 
39  ~HLTTrackWithHits() override {}
40 
41  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
44  desc.add<edm::InputTag>("src", edm::InputTag(""));
45  desc.add<int>("MinN", 0);
46  desc.add<int>("MaxN", 99999);
47  desc.add<int>("MinBPX", 0);
48  desc.add<int>("MinFPX", 0);
49  desc.add<int>("MinPXL", 0);
50  desc.add<double>("MinPT", 0.);
51  descriptions.add("hltTrackWithHits", desc);
52  }
53 
54 private:
56  const edm::EventSetup&,
57  trigger::TriggerFilterObjectWithRefs& filterproduct) const override {
59  iEvent.getByToken(srcToken_, oHandle);
60  int s = oHandle->size();
61  int count = 0;
62  for (int i = 0; i != s; ++i) {
63  const reco::Track& track = (*oHandle)[i];
64  if (track.pt() < MinPT_)
65  continue;
66  const reco::HitPattern& hits = track.hitPattern();
67  if (MinBPX_ > 0 && hits.numberOfValidPixelBarrelHits() >= MinBPX_) {
68  ++count;
69  continue;
70  }
71  if (MinFPX_ > 0 && hits.numberOfValidPixelEndcapHits() >= MinFPX_) {
72  ++count;
73  continue;
74  }
75  if (MinPXL_ > 0 && hits.numberOfValidPixelHits() >= MinPXL_) {
76  ++count;
77  continue;
78  }
79  }
80 
81  bool answer = (count >= minN_ && count <= maxN_);
82  LogDebug("HLTTrackWithHits") << module(iEvent) << " sees: " << s << " objects. Only: " << count
83  << " satisfy the hit requirement. Filter answer is: " << (answer ? "true" : "false")
84  << std::endl;
85  return answer;
86  }
87 
91  double MinPT_;
92 };
93 
94 #endif
bool hltFilter(edm::Event &iEvent, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
int iEvent
Definition: GenABIO.cc:224
int numberOfValidPixelBarrelHits() const
Definition: HitPattern.h:829
double pt() const
track transverse momentum
Definition: TrackBase.h:637
int module(edm::Event const &) const
Definition: HLTFilter.cc:47
HLTTrackWithHits(const edm::ParameterSet &iConfig)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:25
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:504
void add(std::string const &label, ParameterSetDescription const &psetDescription)
int numberOfValidPixelEndcapHits() const
Definition: HitPattern.h:833
edm::EDGetTokenT< reco::TrackCollection > srcToken_
int numberOfValidPixelHits() const
Definition: HitPattern.h:825
~HLTTrackWithHits() override
edm::InputTag src_
#define LogDebug(id)