CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 
31 public:
32 
35 
42  const AnnealingSchedule & ann = GeometricAnnealing(),
43  const LinearizationPointFinder & linP =
45  const VertexUpdator<5> & updator = KalmanVertexUpdator<5>(),
48  const VertexSmoother<5> & smoother = DummyVertexSmoother<5>(),
50 
52 
53  virtual ~AdaptiveVertexFitter();
54 
61  virtual CachingVertex<5> vertex( const std::vector<reco::TransientTrack> & ) const;
62 
71  virtual CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack> & ) const;
72 
76  virtual CachingVertex<5> vertex(const std::vector<RefCountedVertexTrack> &,
77  const reco::BeamSpot & spot ) const;
78 
82  virtual CachingVertex<5> vertex( const std::vector<reco::TransientTrack> &,
83  const GlobalPoint& linPoint ) const;
84 
90  virtual CachingVertex<5> vertex( const std::vector<reco::TransientTrack> &,
91  const GlobalPoint & priorPos,
92  const GlobalError & priorError ) const;
93 
98  virtual CachingVertex<5> vertex(const std::vector<reco::TransientTrack> & tracks,
99  const reco::BeamSpot& beamSpot) const;
100 
101 
106  virtual CachingVertex<5> vertex( const std::vector<RefCountedVertexTrack> &,
107  const GlobalPoint & priorPos,
108  const GlobalError & priorError ) const;
109 
110  AdaptiveVertexFitter * clone() const;
111 
118  void setWeightThreshold ( float w );
119 
131  void setParameters( double maxshift=0.0001, double maxlpshift=0.1,
132  unsigned maxstep=30, double weightthreshold=.001 );
133 
139  void setParameters ( const edm::ParameterSet & );
140 
142 
144 
145 private:
156  std::vector<RefCountedVertexTrack> reLinearizeTracks(
157  const std::vector<RefCountedVertexTrack> & tracks,
158  const CachingVertex<5> & vertex ) const;
159 
160 
165  std::vector<RefCountedVertexTrack> reWeightTracks(
166  const std::vector<RefCountedLinearizedTrackState> &,
167  const CachingVertex<5> & seed ) const;
168 
174  std::vector<RefCountedVertexTrack> reWeightTracks(
175  const std::vector<RefCountedVertexTrack> &,
176  const CachingVertex<5> & seed) const;
177 
178 
184  std::vector<RefCountedVertexTrack> weightTracks(
185  const std::vector<RefCountedLinearizedTrackState> &,
186  const VertexState & seed ) const;
187 
191  std::vector<RefCountedVertexTrack> linearizeTracks(
192  const std::vector<reco::TransientTrack> &,
193  const VertexState & ) const;
197  CachingVertex<5> fit( const std::vector<RefCountedVertexTrack> & tracks,
198  const VertexState & priorSeed,
199  bool withPrior) const;
200 
201  double getWeight ( float chi2 ) const;
202 private:
203  double theMaxShift;
207  mutable int theNr;
208 
216 };
217 
218 #endif
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
CachingVertex< 5 > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState &priorSeed, bool withPrior) const
std::vector< RefCountedVertexTrack > reLinearizeTracks(const std::vector< RefCountedVertexTrack > &tracks, const CachingVertex< 5 > &vertex) const
AdaptiveVertexFitter * clone() const
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) 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())
tuple tracks
Definition: testEve_cfg.py:39
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