CMS 3D CMS Logo

HLTPixelTrackFilter.cc
Go to the documentation of this file.
3 
10 
15 //
16 // class declaration
17 //
18 
20 public:
21  explicit HLTPixelTrackFilter(const edm::ParameterSet&);
22  ~HLTPixelTrackFilter() override;
23  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
24 
25 private:
26  bool filter(edm::Event&, const edm::EventSetup&) override;
27 
28  edm::InputTag inputTag_; // input tag identifying product containing pixel clusters
29  unsigned int min_pixelTracks_; // minimum number of clusters
30  unsigned int max_pixelTracks_; // maximum number of clusters
32 
33 };
34 
35 //
36 // constructors and destructor
37 //
38 
40  inputTag_ (config.getParameter<edm::InputTag>("pixelTracks")),
41  min_pixelTracks_ (config.getParameter<unsigned int>("minPixelTracks")),
42  max_pixelTracks_ (config.getParameter<unsigned int>("maxPixelTracks"))
43 {
44  inputToken_ = consumes< reco::TrackCollection >(inputTag_);
45  LogDebug("") << "Using the " << inputTag_ << " input collection";
46  LogDebug("") << "Requesting at least " << min_pixelTracks_ << " PixelTracks";
47  if(max_pixelTracks_ > 0)
48  LogDebug("") << "...but no more than " << max_pixelTracks_ << " PixelTracks";
49 }
50 
52 
53 void
56  desc.add<edm::InputTag>("pixelTracks",edm::InputTag("hltPixelTracks"));
57  desc.add<unsigned int>("minPixelTracks",0);
58  desc.add<unsigned int>("maxPixelTracks",0);
59  descriptions.add("hltPixelTrackFilter",desc);
60 
61 }
62 
63 //
64 // member functions
65 //
66 // ------------ method called to produce the data ------------
68 {
69  // get hold of products from Event
71  event.getByToken(inputToken_, trackColl);
72 
73  unsigned int numTracks = trackColl->size();
74  LogDebug("") << "Number of tracks accepted: " << numTracks;
75  bool accept = (numTracks >= min_pixelTracks_);
76  if(max_pixelTracks_ > 0)
77  accept &= (numTracks <= max_pixelTracks_);
78  // return with final filter decision
79  return accept;
80 }
81 
82 // define as a framework module
#define LogDebug(id)
HLTPixelTrackFilter(const edm::ParameterSet &)
Definition: config.py:1
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool filter(edm::Event &, const edm::EventSetup &) override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
~HLTPixelTrackFilter() override
Definition: event.py:1
edm::EDGetTokenT< reco::TrackCollection > inputToken_