File indexing completed on 2024-04-06 12:29:09
0001 #ifndef KalmanVertexFitter_H
0002 #define KalmanVertexFitter_H
0003
0004 #include "RecoVertex/VertexTools/interface/SequentialVertexFitter.h"
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022 class KalmanVertexFitter : public VertexFitter<5> {
0023 public:
0024
0025
0026
0027
0028
0029
0030 KalmanVertexFitter(bool useSmoothing = false);
0031
0032
0033
0034
0035
0036 KalmanVertexFitter(const edm::ParameterSet& pSet, bool useSmoothing = false);
0037
0038 KalmanVertexFitter(const KalmanVertexFitter& other) : theSequentialFitter(other.theSequentialFitter->clone()) {}
0039
0040 ~KalmanVertexFitter() override { delete theSequentialFitter; }
0041
0042 KalmanVertexFitter* clone() const override { return new KalmanVertexFitter(*this); }
0043
0044 public:
0045 typedef CachingVertex<5>::RefCountedVertexTrack RefCountedVertexTrack;
0046
0047
0048
0049 inline CachingVertex<5> vertex(const std::vector<reco::TransientTrack>& tracks) const override {
0050 return theSequentialFitter->vertex(tracks);
0051 }
0052
0053
0054
0055 inline CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack>& tracks) const override {
0056 return theSequentialFitter->vertex(tracks);
0057 }
0058
0059 inline CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack>& tracks,
0060 const reco::BeamSpot& spot) const override {
0061 return theSequentialFitter->vertex(tracks, spot);
0062 }
0063
0064
0065
0066
0067 inline CachingVertex<5> vertex(const std::vector<reco::TransientTrack>& tracks,
0068 const GlobalPoint& linPoint) const override {
0069 return theSequentialFitter->vertex(tracks, linPoint);
0070 }
0071
0072
0073
0074
0075
0076
0077 inline CachingVertex<5> vertex(const std::vector<reco::TransientTrack>& tracks,
0078 const GlobalPoint& priorPos,
0079 const GlobalError& priorError) const override {
0080 return theSequentialFitter->vertex(tracks, priorPos, priorError);
0081 }
0082
0083
0084
0085
0086
0087 inline CachingVertex<5> vertex(const std::vector<reco::TransientTrack>& tracks,
0088 const reco::BeamSpot& beamSpot) const override {
0089 return theSequentialFitter->vertex(tracks, beamSpot);
0090 }
0091
0092
0093
0094
0095
0096 inline CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack>& tracks,
0097 const GlobalPoint& priorPos,
0098 const GlobalError& priorError) const override {
0099 return theSequentialFitter->vertex(tracks, priorPos, priorError);
0100 }
0101
0102
0103
0104
0105
0106 private:
0107 void setup(const edm::ParameterSet& pSet, bool useSmoothing);
0108
0109 edm::ParameterSet defaultParameters() const;
0110
0111 const SequentialVertexFitter<5>* theSequentialFitter;
0112 };
0113
0114 #endif