CMS 3D CMS Logo

AdaptiveGsfVertexFitter.h
Go to the documentation of this file.
1 #ifndef AdaptiveGsfVertexFitter_H
2 #define AdaptiveGsfVertexFitter_H
3 
8 
17 
18 public:
19 
21 
29 
30  ~AdaptiveGsfVertexFitter() override;
31 
37 
38  AdaptiveGsfVertexFitter * clone() const override {
39  return new AdaptiveGsfVertexFitter(* this);
40  }
41 
42 public:
43 
46  inline CachingVertex<5>
47  vertex(const std::vector<reco::TransientTrack> & tracks) const override
48  {
49  return theFitter->vertex(tracks);
50  }
51 
54  inline CachingVertex<5>
55  vertex(const std::vector<RefCountedVertexTrack> & tracks) const override
56  {
57  return theFitter->vertex(tracks);
58  }
59 
63  inline CachingVertex<5>
64  vertex(const std::vector<reco::TransientTrack> & tracks,
65  const GlobalPoint& linPoint) const override
66  {
67  return theFitter->vertex(tracks, linPoint);
68  }
69 
74  inline CachingVertex<5>
75  vertex(const std::vector<reco::TransientTrack> & tracks, const reco::BeamSpot& beamSpot) const override
76  {
77  return theFitter->vertex(tracks, beamSpot);
78  }
79 
80 
86  inline CachingVertex<5>
87  vertex(const std::vector<reco::TransientTrack> & tracks,
88  const GlobalPoint& priorPos,
89  const GlobalError& priorError) const override
90  {
91  return theFitter->vertex(tracks, priorPos, priorError);
92  }
93 
94  inline CachingVertex<5>
95  vertex(const std::vector<RefCountedVertexTrack> & tracks,
96  const reco::BeamSpot & spot ) const override
97  {
98  return theFitter->vertex(tracks, spot );
99  }
100 
105  inline CachingVertex<5>
106  vertex(const std::vector<RefCountedVertexTrack> & tracks,
107  const GlobalPoint& priorPos,
108  const GlobalError& priorError) const override
109  {
110  return theFitter->vertex(tracks, priorPos, priorError);
111  }
112 
113 private:
114 
116 };
117 
118 #endif
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks, const GlobalPoint &linPoint) const override
AdaptiveVertexFitter * theFitter
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const override
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks, const GlobalPoint &priorPos, const GlobalError &priorError) const override
CachingVertex< 5 > vertex(const std::vector< RefCountedVertexTrack > &tracks, const GlobalPoint &priorPos, const GlobalError &priorError) const override
CachingVertex< 5 > vertex(const std::vector< RefCountedVertexTrack > &tracks, const reco::BeamSpot &spot) const override
AdaptiveGsfVertexFitter(const edm::ParameterSet &pSet, const LinearizationPointFinder &linP=DefaultLinearizationPointFinder())
CachingVertex< 5 > vertex(const std::vector< RefCountedVertexTrack > &tracks) const override
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
AdaptiveGsfVertexFitter * clone() const override
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks, const reco::BeamSpot &beamSpot) const override
CachingVertex< 5 >::RefCountedVertexTrack RefCountedVertexTrack