Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 14:20:42

0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/stream/EDFilter.h"
0003 
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "DataFormats/Common/interface/Handle.h"
0008 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0009 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0010 
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014 #include "DataFormats/TrackReco/interface/Track.h"
0015 //
0016 // class declaration
0017 //
0018 
0019 class HLTPixelTrackFilter : public edm::stream::EDFilter<> {
0020 public:
0021   explicit HLTPixelTrackFilter(const edm::ParameterSet&);
0022   ~HLTPixelTrackFilter() override;
0023   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0024 
0025 private:
0026   bool filter(edm::Event&, const edm::EventSetup&) override;
0027 
0028   edm::InputTag inputTag_;        // input tag identifying product containing pixel clusters
0029   unsigned int min_pixelTracks_;  // minimum number of clusters
0030   unsigned int max_pixelTracks_;  // maximum number of clusters
0031   edm::EDGetTokenT<reco::TrackCollection> inputToken_;
0032 };
0033 
0034 //
0035 // constructors and destructor
0036 //
0037 
0038 HLTPixelTrackFilter::HLTPixelTrackFilter(const edm::ParameterSet& config)
0039     : inputTag_(config.getParameter<edm::InputTag>("pixelTracks")),
0040       min_pixelTracks_(config.getParameter<unsigned int>("minPixelTracks")),
0041       max_pixelTracks_(config.getParameter<unsigned int>("maxPixelTracks")) {
0042   inputToken_ = consumes<reco::TrackCollection>(inputTag_);
0043   LogDebug("") << "Using the " << inputTag_ << " input collection";
0044   LogDebug("") << "Requesting at least " << min_pixelTracks_ << " PixelTracks";
0045   if (max_pixelTracks_ > 0)
0046     LogDebug("") << "...but no more than " << max_pixelTracks_ << " PixelTracks";
0047 }
0048 
0049 HLTPixelTrackFilter::~HLTPixelTrackFilter() = default;
0050 
0051 void HLTPixelTrackFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0052   edm::ParameterSetDescription desc;
0053   desc.add<edm::InputTag>("pixelTracks", edm::InputTag("hltPixelTracks"));
0054   desc.add<unsigned int>("minPixelTracks", 0);
0055   desc.add<unsigned int>("maxPixelTracks", 0);
0056   descriptions.add("hltPixelTrackFilter", desc);
0057 }
0058 
0059 //
0060 // member functions
0061 //
0062 // ------------ method called to produce the data  ------------
0063 bool HLTPixelTrackFilter::filter(edm::Event& event, const edm::EventSetup& iSetup) {
0064   // get hold of products from Event
0065   edm::Handle<reco::TrackCollection> trackColl;
0066   event.getByToken(inputToken_, trackColl);
0067 
0068   unsigned int numTracks = trackColl->size();
0069   LogDebug("") << "Number of tracks accepted: " << numTracks;
0070   bool accept = (numTracks >= min_pixelTracks_);
0071   if (max_pixelTracks_ > 0)
0072     accept &= (numTracks <= max_pixelTracks_);
0073   // return with final filter decision
0074   return accept;
0075 }
0076 
0077 // define as a framework module
0078 #include "FWCore/Framework/interface/MakerMacros.h"
0079 DEFINE_FWK_MODULE(HLTPixelTrackFilter);