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
13 
16 
18 
23 
24 class HLTTrackWithHits : public HLTFilter {
25 public:
26  explicit HLTTrackWithHits(const edm::ParameterSet& iConfig)
27  : HLTFilter(iConfig),
28  src_(iConfig.getParameter<edm::InputTag>("src")),
29  minN_(iConfig.getParameter<int>("MinN")),
30  maxN_(iConfig.getParameter<int>("MaxN")),
31  MinBPX_(iConfig.getParameter<int>("MinBPX")),
32  MinFPX_(iConfig.getParameter<int>("MinFPX")),
33  MinPXL_(iConfig.getParameter<int>("MinPXL")),
34  MinPT_(iConfig.getParameter<double>("MinPT")) {
35  srcToken_ = consumes<reco::TrackCollection>(src_);
36  }
37 
38  ~HLTTrackWithHits() override {}
39 
40  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
43  desc.add<edm::InputTag>("src", edm::InputTag(""));
44  desc.add<int>("MinN", 0);
45  desc.add<int>("MaxN", 99999);
46  desc.add<int>("MinBPX", 0);
47  desc.add<int>("MinFPX", 0);
48  desc.add<int>("MinPXL", 0);
49  desc.add<double>("MinPT", 0.);
50  descriptions.add("hltTrackWithHits", desc);
51  }
52 
53 private:
55  const edm::EventSetup&,
56  trigger::TriggerFilterObjectWithRefs& filterproduct) const override {
58  iEvent.getByToken(srcToken_, oHandle);
59  int s = oHandle->size();
60  int count = 0;
61  for (int i = 0; i != s; ++i) {
62  const reco::Track& track = (*oHandle)[i];
63  if (track.pt() < MinPT_)
64  continue;
65  const reco::HitPattern& hits = track.hitPattern();
66  if (MinBPX_ > 0 && hits.numberOfValidPixelBarrelHits() >= MinBPX_) {
67  ++count;
68  continue;
69  }
70  if (MinFPX_ > 0 && hits.numberOfValidPixelEndcapHits() >= MinFPX_) {
71  ++count;
72  continue;
73  }
74  if (MinPXL_ > 0 && hits.numberOfValidPixelHits() >= MinPXL_) {
75  ++count;
76  continue;
77  }
78  }
79 
80  bool answer = (count >= minN_ && count <= maxN_);
81  LogDebug("HLTTrackWithHits") << module(iEvent) << " sees: " << s << " objects. Only: " << count
82  << " satisfy the hit requirement. Filter answer is: " << (answer ? "true" : "false")
83  << std::endl;
84  return answer;
85  }
86 
90  double MinPT_;
91 };
92 
93 #endif
bool hltFilter(edm::Event &iEvent, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int module(edm::Event const &) const
Definition: HLTFilter.cc:47
int iEvent
Definition: GenABIO.cc:224
HLTTrackWithHits(const edm::ParameterSet &iConfig)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:25
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
edm::EDGetTokenT< reco::TrackCollection > srcToken_
~HLTTrackWithHits() override
edm::InputTag src_
#define LogDebug(id)