File indexing completed on 2024-07-10 02:34:48
0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "HLTrigger/HLTcore/interface/HLTFilter.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/RecoCandidate/interface/RecoChargedCandidate.h"
0015
0016
0017
0018
0019 class HLTPixelTrackFilter : public HLTFilter {
0020 public:
0021 explicit HLTPixelTrackFilter(const edm::ParameterSet&);
0022 ~HLTPixelTrackFilter() override;
0023 bool hltFilter(edm::Event&,
0024 const edm::EventSetup&,
0025 trigger::TriggerFilterObjectWithRefs& filterproduct) const override;
0026 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0027
0028 private:
0029 edm::InputTag inputTag_;
0030 unsigned int min_pixelTracks_;
0031 unsigned int max_pixelTracks_;
0032 edm::EDGetTokenT<reco::RecoChargedCandidateCollection> inputToken_;
0033 };
0034
0035
0036
0037
0038
0039 HLTPixelTrackFilter::HLTPixelTrackFilter(const edm::ParameterSet& config)
0040 : HLTFilter(config),
0041 inputTag_(config.getParameter<edm::InputTag>("pixelTracks")),
0042 min_pixelTracks_(config.getParameter<unsigned int>("minPixelTracks")),
0043 max_pixelTracks_(config.getParameter<unsigned int>("maxPixelTracks")) {
0044 inputToken_ = consumes<reco::RecoChargedCandidateCollection>(inputTag_);
0045 LogDebug("") << "Using the " << inputTag_ << " input collection";
0046 LogDebug("") << "Requesting at least " << min_pixelTracks_ << " PixelTracks";
0047 if (max_pixelTracks_ > 0)
0048 LogDebug("") << "...but no more than " << max_pixelTracks_ << " PixelTracks";
0049 }
0050
0051 HLTPixelTrackFilter::~HLTPixelTrackFilter() = default;
0052
0053 void HLTPixelTrackFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0054 edm::ParameterSetDescription desc;
0055 makeHLTFilterDescription(desc);
0056 desc.add<edm::InputTag>("pixelTracks", edm::InputTag("hltPixelTracks"));
0057 desc.add<unsigned int>("minPixelTracks", 0);
0058 desc.add<unsigned int>("maxPixelTracks", 0);
0059 descriptions.add("hltPixelTrackFilter", desc);
0060 }
0061
0062
0063
0064
0065
0066 bool HLTPixelTrackFilter::hltFilter(edm::Event& event,
0067 const edm::EventSetup& iSetup,
0068 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0069
0070
0071
0072 if (saveTags())
0073 filterproduct.addCollectionTag(inputTag_);
0074
0075
0076 const auto& trackColl = event.getHandle(inputToken_);
0077
0078 unsigned int numTracks = trackColl->size();
0079 for (size_t i = 0; i < numTracks; i++)
0080 filterproduct.addObject(trigger::TriggerTrack, reco::RecoChargedCandidateRef(trackColl, i));
0081
0082 LogDebug("") << "Number of tracks accepted: " << numTracks;
0083 bool accept = (numTracks >= min_pixelTracks_);
0084 if (max_pixelTracks_ > 0)
0085 accept &= (numTracks <= max_pixelTracks_);
0086
0087 return accept;
0088 }
0089
0090
0091 #include "FWCore/Framework/interface/MakerMacros.h"
0092 DEFINE_FWK_MODULE(HLTPixelTrackFilter);