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)
30  : HLTFilter(iConfig),
31  src_(iConfig.getParameter<edm::InputTag>("src")),
32  minN_(iConfig.getParameter<int>("MinN")),
33  maxN_(iConfig.getParameter<int>("MaxN")),
34  MinBPX_(iConfig.getParameter<int>("MinBPX")),
35  MinFPX_(iConfig.getParameter<int>("MinFPX")),
36  MinPXL_(iConfig.getParameter<int>("MinPXL")),
37  MinPT_(iConfig.getParameter<double>("MinPT")) {
38  srcToken_ = consumes<reco::TrackCollection>(src_);
39  }
40 
41  ~HLTTrackWithHits() override {}
42 
43  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
46  desc.add<edm::InputTag>("src", edm::InputTag(""));
47  desc.add<int>("MinN", 0);
48  desc.add<int>("MaxN", 99999);
49  desc.add<int>("MinBPX", 0);
50  desc.add<int>("MinFPX", 0);
51  desc.add<int>("MinPXL", 0);
52  desc.add<double>("MinPT", 0.);
53  descriptions.add("hltTrackWithHits", desc);
54  }
55 
56 private:
58  const edm::EventSetup&,
59  trigger::TriggerFilterObjectWithRefs& filterproduct) const override {
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_)
67  continue;
68  const reco::HitPattern& hits = track.hitPattern();
69  if (MinBPX_ > 0 && hits.numberOfValidPixelBarrelHits() >= MinBPX_) {
70  ++count;
71  continue;
72  }
73  if (MinFPX_ > 0 && hits.numberOfValidPixelEndcapHits() >= MinFPX_) {
74  ++count;
75  continue;
76  }
77  if (MinPXL_ > 0 && hits.numberOfValidPixelHits() >= MinPXL_) {
78  ++count;
79  continue;
80  }
81  }
82 
83  bool answer = (count >= minN_ && count <= maxN_);
84  LogDebug("HLTTrackWithHits") << module(iEvent) << " sees: " << s << " objects. Only: " << count
85  << " satisfy the hit requirement. Filter answer is: " << (answer ? "true" : "false")
86  << std::endl;
87  return answer;
88  }
89 
93  double MinPT_;
94 };
95 
96 #endif
#define LogDebug(id)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
int iEvent
Definition: GenABIO.cc:224
int numberOfValidPixelBarrelHits() const
Definition: HitPattern.h:805
double pt() const
track transverse momentum
Definition: TrackBase.h:602
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:483
void add(std::string const &label, ParameterSetDescription const &psetDescription)
int numberOfValidPixelEndcapHits() const
Definition: HitPattern.h:809
HLT enums.
edm::EDGetTokenT< reco::TrackCollection > srcToken_
int numberOfValidPixelHits() const
Definition: HitPattern.h:801
bool hltFilter(edm::Event &iEvent, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
~HLTTrackWithHits() override
edm::InputTag src_
answer
Definition: submit.py:45