File indexing completed on 2023-03-17 11:23:28
0001 #ifndef HITrackFilterForPVFinding_h
0002 #define HITrackFilterForPVFinding_h
0003
0004
0005
0006
0007
0008
0009
0010
0011 #include "RecoVertex/PrimaryVertexProducer/interface/TrackFilterForPVFinding.h"
0012
0013 class HITrackFilterForPVFinding : public TrackFilterForPVFinding {
0014 private:
0015 unsigned int NumTracksThreshold_;
0016 unsigned int MaxNumTracksThreshold_;
0017 double minPtTight_;
0018
0019 public:
0020 HITrackFilterForPVFinding(const edm::ParameterSet& conf) : TrackFilterForPVFinding(conf) {
0021 NumTracksThreshold_ = conf.getParameter<int>("numTracksThreshold");
0022 MaxNumTracksThreshold_ = conf.getParameter<int>("maxNumTracksThreshold");
0023 minPtTight_ = conf.getParameter<double>("minPtTight");
0024
0025 }
0026
0027
0028 std::vector<reco::TransientTrack> select(const std::vector<reco::TransientTrack>& tracks) const override {
0029 std::vector<reco::TransientTrack> seltks = TrackFilterForPVFinding::select(tracks);
0030 if (seltks.size() < NumTracksThreshold_) {
0031 return tracks;
0032 } else if (seltks.size() > MaxNumTracksThreshold_) {
0033 std::vector<reco::TransientTrack> seltksTight = TrackFilterForPVFinding::selectTight(tracks, minPtTight_);
0034 if (seltksTight.size() >= NumTracksThreshold_)
0035 return seltksTight;
0036 }
0037
0038 return seltks;
0039 }
0040
0041 static void fillPSetDescription(edm::ParameterSetDescription& desc) {
0042 TrackFilterForPVFinding::fillPSetDescription(desc);
0043 desc.add<int>("numTracksThreshold", 0);
0044 desc.add<int>("maxNumTracksThreshold", std::numeric_limits<int>::max());
0045 desc.add<double>("minPtTight", 0.0);
0046 }
0047 };
0048
0049 #endif