File indexing completed on 2024-04-06 12:29:18
0001 #ifndef _ConfigurableTrimmedVertexFinder_H_
0002 #define _ConfigurableTrimmedVertexFinder_H_
0003
0004 #include "RecoVertex/VertexPrimitives/interface/VertexReconstructor.h"
0005 #include "RecoVertex/TrimmedKalmanVertexFinder/interface/TrimmedVertexFinder.h"
0006 #include "RecoVertex/TrimmedKalmanVertexFinder/interface/TrimmedTrackFilter.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include <vector>
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044 class ConfigurableTrimmedVertexFinder : public VertexReconstructor {
0045 public:
0046 ConfigurableTrimmedVertexFinder(const VertexFitter<5>* vf,
0047 const VertexUpdator<5>* vu,
0048 const VertexTrackCompatibilityEstimator<5>* ve);
0049
0050 ~ConfigurableTrimmedVertexFinder() override {}
0051
0052 std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack>& tracks) const override;
0053
0054 std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack>& tracks,
0055 const reco::BeamSpot& spot) const override;
0056
0057 std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack>& tracks,
0058 std::vector<reco::TransientTrack>& unused,
0059 const reco::BeamSpot& spot,
0060 bool use_spot) const;
0061
0062
0063
0064 float ptCut() const { return theFilter.ptCut(); }
0065 const TrimmedTrackFilter& trackFilter() const { return theFilter; }
0066 float trackCompatibilityCut() const { return theTrackCompatibilityToPV; }
0067 float trackCompatibilityToSV() const { return theTrackCompatibilityToSV; }
0068 float vertexFitProbabilityCut() const { return theVtxFitProbCut; }
0069 int maxNbOfVertices() const { return theMaxNbOfVertices; }
0070
0071
0072
0073 void setParameters(const edm::ParameterSet&);
0074
0075 void setPtCut(float cut) { theFilter.setPtCut(cut); }
0076 void setTrackCompatibilityCut(float cut) { theTrackCompatibilityToPV = cut; }
0077 void setTrackCompatibilityToSV(float cut) { theTrackCompatibilityToSV = cut; }
0078 void setVertexFitProbabilityCut(float cut) { theVtxFitProbCut = cut; }
0079 void setMaxNbOfVertices(int max) { theMaxNbOfVertices = max; }
0080
0081
0082
0083 ConfigurableTrimmedVertexFinder* clone() const override { return new ConfigurableTrimmedVertexFinder(*this); }
0084
0085 protected:
0086 virtual void resetEvent(const std::vector<reco::TransientTrack>& tracks) const {}
0087
0088 virtual void analyseInputTracks(const std::vector<reco::TransientTrack>& tracks) const {}
0089
0090 virtual void analyseClusterFinder(const std::vector<TransientVertex>& vts,
0091 const std::vector<reco::TransientTrack>& remain) const {}
0092
0093 virtual void analyseVertexCandidates(const std::vector<TransientVertex>& vts) const {}
0094
0095 virtual void analyseFoundVertices(const std::vector<TransientVertex>& vts) const {}
0096
0097 private:
0098 using VertexReconstructor::vertices;
0099
0100
0101
0102
0103
0104
0105 std::vector<TransientVertex> vertexCandidates(const std::vector<reco::TransientTrack>& tracks,
0106 std::vector<reco::TransientTrack>& unused,
0107 const reco::BeamSpot& spot,
0108 bool use_spot) const;
0109
0110
0111 std::vector<TransientVertex> clean(const std::vector<TransientVertex>& candidates) const;
0112
0113
0114 mutable TrimmedVertexFinder theClusterFinder;
0115 float theVtxFitProbCut;
0116 float theTrackCompatibilityToPV;
0117 float theTrackCompatibilityToSV;
0118 int theMaxNbOfVertices;
0119 TrimmedTrackFilter theFilter;
0120 };
0121
0122 #endif