File indexing completed on 2024-04-06 12:18:43
0001 #ifndef HLTrigger_HLTTrackWithHits_H
0002
0003
0004
0005
0006
0007
0008
0009 #include <memory>
0010
0011
0012 #include "FWCore/Framework/interface/Frameworkfwd.h"
0013
0014 #include "FWCore/Framework/interface/Event.h"
0015 #include "FWCore/Framework/interface/MakerMacros.h"
0016
0017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0018
0019 #include "HLTrigger/HLTcore/interface/HLTFilter.h"
0020 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0021 #include "DataFormats/TrackReco/interface/Track.h"
0022 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0023
0024 class HLTTrackWithHits : public HLTFilter {
0025 public:
0026 explicit HLTTrackWithHits(const edm::ParameterSet& iConfig)
0027 : HLTFilter(iConfig),
0028 src_(iConfig.getParameter<edm::InputTag>("src")),
0029 minN_(iConfig.getParameter<int>("MinN")),
0030 maxN_(iConfig.getParameter<int>("MaxN")),
0031 MinBPX_(iConfig.getParameter<int>("MinBPX")),
0032 MinFPX_(iConfig.getParameter<int>("MinFPX")),
0033 MinPXL_(iConfig.getParameter<int>("MinPXL")),
0034 MinPT_(iConfig.getParameter<double>("MinPT")) {
0035 srcToken_ = consumes<reco::TrackCollection>(src_);
0036 }
0037
0038 ~HLTTrackWithHits() override {}
0039
0040 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0041 edm::ParameterSetDescription desc;
0042 makeHLTFilterDescription(desc);
0043 desc.add<edm::InputTag>("src", edm::InputTag(""));
0044 desc.add<int>("MinN", 0);
0045 desc.add<int>("MaxN", 99999);
0046 desc.add<int>("MinBPX", 0);
0047 desc.add<int>("MinFPX", 0);
0048 desc.add<int>("MinPXL", 0);
0049 desc.add<double>("MinPT", 0.);
0050 descriptions.add("hltTrackWithHits", desc);
0051 }
0052
0053 private:
0054 bool hltFilter(edm::Event& iEvent,
0055 const edm::EventSetup&,
0056 trigger::TriggerFilterObjectWithRefs& filterproduct) const override {
0057 edm::Handle<reco::TrackCollection> oHandle;
0058 iEvent.getByToken(srcToken_, oHandle);
0059 int s = oHandle->size();
0060 int count = 0;
0061 for (int i = 0; i != s; ++i) {
0062 const reco::Track& track = (*oHandle)[i];
0063 if (track.pt() < MinPT_)
0064 continue;
0065 const reco::HitPattern& hits = track.hitPattern();
0066 if (MinBPX_ > 0 && hits.numberOfValidPixelBarrelHits() >= MinBPX_) {
0067 ++count;
0068 continue;
0069 }
0070 if (MinFPX_ > 0 && hits.numberOfValidPixelEndcapHits() >= MinFPX_) {
0071 ++count;
0072 continue;
0073 }
0074 if (MinPXL_ > 0 && hits.numberOfValidPixelHits() >= MinPXL_) {
0075 ++count;
0076 continue;
0077 }
0078 }
0079
0080 bool answer = (count >= minN_ && count <= maxN_);
0081 LogDebug("HLTTrackWithHits") << module(iEvent) << " sees: " << s << " objects. Only: " << count
0082 << " satisfy the hit requirement. Filter answer is: " << (answer ? "true" : "false")
0083 << std::endl;
0084 return answer;
0085 }
0086
0087 edm::InputTag src_;
0088 edm::EDGetTokenT<reco::TrackCollection> srcToken_;
0089 int minN_, maxN_, MinBPX_, MinFPX_, MinPXL_;
0090 double MinPT_;
0091 };
0092
0093 #endif