CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
AdaptiveVertexFitter.h
Go to the documentation of this file.
1 #ifndef AdaptiveVertexFitter_H
2 #define AdaptiveVertexFitter_H
3 
9 // #include "RecoVertex/AdaptiveVertexFit/interface/KalmanVertexSmoother.h"
15 
30 public:
33 
43  const VertexSmoother<5> &smoother = DummyVertexSmoother<5>(),
45 
47 
48  ~AdaptiveVertexFitter() override;
49 
56  CachingVertex<5> vertex(const std::vector<reco::TransientTrack> &) const override;
57 
66  CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack> &) const override;
67 
71  CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack> &, const reco::BeamSpot &spot) const override;
72 
76  CachingVertex<5> vertex(const std::vector<reco::TransientTrack> &, const GlobalPoint &linPoint) const override;
77 
83  CachingVertex<5> vertex(const std::vector<reco::TransientTrack> &,
84  const GlobalPoint &priorPos,
85  const GlobalError &priorError) const override;
86 
91  CachingVertex<5> vertex(const std::vector<reco::TransientTrack> &tracks,
92  const reco::BeamSpot &beamSpot) const override;
93 
98  CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack> &,
99  const GlobalPoint &priorPos,
100  const GlobalError &priorError) const override;
101 
102  AdaptiveVertexFitter *clone() const override;
103 
110  void setWeightThreshold(float w);
111 
123  void setParameters(double maxshift = 0.0001,
124  double maxlpshift = 0.1,
125  unsigned maxstep = 30,
126  double weightthreshold = .001);
127 
133  void setParameters(const edm::ParameterSet &);
134 
136 
138 
139 private:
150  std::vector<RefCountedVertexTrack> reLinearizeTracks(const std::vector<RefCountedVertexTrack> &tracks,
151  const CachingVertex<5> &vertex) const;
152 
157  std::vector<RefCountedVertexTrack> reWeightTracks(const std::vector<RefCountedLinearizedTrackState> &,
158  const CachingVertex<5> &seed) const;
159 
165  std::vector<RefCountedVertexTrack> reWeightTracks(const std::vector<RefCountedVertexTrack> &,
166  const CachingVertex<5> &seed) const;
167 
173  std::vector<RefCountedVertexTrack> weightTracks(const std::vector<RefCountedLinearizedTrackState> &,
174  const VertexState &seed) const;
175 
179  std::vector<RefCountedVertexTrack> linearizeTracks(const std::vector<reco::TransientTrack> &,
180  const VertexState &) const;
184  CachingVertex<5> fit(const std::vector<RefCountedVertexTrack> &tracks,
185  const VertexState &priorSeed,
186  bool withPrior) const;
187 
188  double getWeight(float chi2) const;
189 
190 private:
191  double theMaxShift;
195  mutable int theNr;
196 
204 };
205 
206 #endif
AdaptiveVertexFitter * clone() const override
std::vector< RefCountedVertexTrack > reWeightTracks(const std::vector< RefCountedLinearizedTrackState > &, const CachingVertex< 5 > &seed) const
const double w
Definition: UKUtility.cc:23
LinearizationPointFinder * theLinP
VertexUpdator< 5 > * theUpdator
list original
Definition: definitions.py:57
auto const & tracks
cannot be loose
CachingVertex< 5 > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState &priorSeed, bool withPrior) const
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const override
std::vector< RefCountedVertexTrack > reLinearizeTracks(const std::vector< RefCountedVertexTrack > &tracks, const CachingVertex< 5 > &vertex) const
VertexSmoother< 5 > * theSmoother
ReferenceCountingPointer< VertexTrack< 5 > > RefCountedVertexTrack
AnnealingSchedule * theAssProbComputer
void setParameters(double maxshift=0.0001, double maxlpshift=0.1, unsigned maxstep=30, double weightthreshold=.001)
double getWeight(float chi2) const
AdaptiveVertexFitter(const AnnealingSchedule &ann=GeometricAnnealing(), const LinearizationPointFinder &linP=DefaultLinearizationPointFinder(), const VertexUpdator< 5 > &updator=KalmanVertexUpdator< 5 >(), const VertexTrackCompatibilityEstimator< 5 > &estor=KalmanVertexTrackCompatibilityEstimator< 5 >(), const VertexSmoother< 5 > &smoother=DummyVertexSmoother< 5 >(), const AbstractLTSFactory< 5 > &ltsf=LinearizedTrackStateFactory())
bool gsfIntermediarySmoothing() const
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &, const VertexState &) const
ReferenceCountingPointer< LinearizedTrackState< 5 > > RefCountedLinearizedTrackState
std::vector< RefCountedVertexTrack > weightTracks(const std::vector< RefCountedLinearizedTrackState > &, const VertexState &seed) const
void gsfIntermediarySmoothing(bool sm)
const AbstractLTSFactory< 5 > * theLinTrkFactory
VertexTrackCompatibilityEstimator< 5 > * theComp