Go to the documentation of this file.00001 #ifndef SequentialVertexFitter_H
00002 #define SequentialVertexFitter_H
00003
00004 #include "RecoVertex/VertexPrimitives/interface/VertexFitter.h"
00005 #include "RecoVertex/VertexTools/interface/LinearizationPointFinder.h"
00006 #include "RecoVertex/VertexPrimitives/interface/VertexUpdator.h"
00007 #include "RecoVertex/VertexPrimitives/interface/VertexSmoother.h"
00008 #include "RecoVertex/VertexTools/interface/LinearizedTrackStateFactory.h"
00009 #include "RecoVertex/VertexTools/interface/VertexTrackFactory.h"
00010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00011 #include <cmath>
00012
00013
00014
00033 template <unsigned int N>
00034 class SequentialVertexFitter : public VertexFitter<N> {
00035
00036 public:
00037
00038 typedef ReferenceCountingPointer<RefittedTrackState<N> > RefCountedRefittedTrackState;
00039 typedef ReferenceCountingPointer<VertexTrack<N> > RefCountedVertexTrack;
00040 typedef ReferenceCountingPointer<LinearizedTrackState<N> > RefCountedLinearizedTrackState;
00041
00042
00049 SequentialVertexFitter(const LinearizationPointFinder & linP,
00050 const VertexUpdator<N> & updator, const VertexSmoother<N> & smoother,
00051 const AbstractLTSFactory<N> & ltsf);
00052
00057 SequentialVertexFitter(const edm::ParameterSet& pSet,
00058 const LinearizationPointFinder & linP,
00059 const VertexUpdator<N> & updator, const VertexSmoother<N> & smoother,
00060 const AbstractLTSFactory<N> & ltsf);
00061
00066 SequentialVertexFitter(const SequentialVertexFitter & original);
00067
00068
00069 virtual ~SequentialVertexFitter();
00070
00071
00078 void setMaximumDistance(float maxShift) {theMaxShift = maxShift;}
00079
00080
00085 void setMaximumNumberOfIterations(int maxIterations)
00086 {theMaxStep = maxIterations;}
00087
00095 virtual CachingVertex<N> vertex(const std::vector<reco::TransientTrack> & tracks) const;
00096
00106 virtual CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack> & tracks) const;
00107
00111 virtual CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack> & tracks, const reco::BeamSpot & spot ) const;
00112
00113
00116 virtual CachingVertex<N> vertex(const std::vector<reco::TransientTrack> & tracks,
00117 const GlobalPoint& linPoint) const;
00118
00123 virtual CachingVertex<N> vertex(const std::vector<reco::TransientTrack> & tracks,
00124 const reco::BeamSpot& beamSpot) const;
00125
00126
00132 virtual CachingVertex<N> vertex(const std::vector<reco::TransientTrack> & tracks,
00133 const GlobalPoint& priorPos,
00134 const GlobalError& priorError) const;
00135
00140 virtual CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack> & tracks,
00141 const GlobalPoint& priorPos,
00142 const GlobalError& priorError) const;
00143
00144
00145
00156
00157
00161 const LinearizationPointFinder * linearizationPointFinder() const
00162 {return theLinP;}
00163
00164 const VertexUpdator<N> * vertexUpdator() const
00165 {return theUpdator;}
00166
00167 const VertexSmoother<N> * vertexSmoother() const
00168 {return theSmoother;}
00169
00170 const float maxShift() const
00171 {return theMaxShift;}
00172
00173 const int maxStep() const
00174 {return theMaxStep;}
00175
00176 const edm::ParameterSet parameterSet() const
00177 {return thePSet;}
00178
00179 SequentialVertexFitter * clone() const {
00180 return new SequentialVertexFitter(* this);
00181 }
00182
00183 const AbstractLTSFactory<N> * linearizedTrackStateFactory() const
00184 { return theLTrackFactory;}
00185
00186 protected:
00187
00192 SequentialVertexFitter() {}
00193
00194 private:
00195
00205 CachingVertex<N> fit(const std::vector<RefCountedVertexTrack> & tracks,
00206 const VertexState priorVertex, bool withPrior) const;
00207
00215 std::vector<RefCountedVertexTrack> linearizeTracks(const std::vector<reco::TransientTrack> & tracks,
00216 const VertexState state) const;
00217
00228 std::vector<RefCountedVertexTrack> reLinearizeTracks(
00229 const std::vector<RefCountedVertexTrack> & tracks,
00230 const VertexState state) const;
00231
00232
00237 void readParameters();
00238 void setDefaultParameters();
00239
00243 inline bool hasNan(const GlobalPoint& point) const {
00244 using namespace std;
00245 return (isnan(point.x())|| isnan(point.y()) || isnan(point.z()));
00246 }
00247
00248
00249 float theMaxShift;
00250 int theMaxStep;
00251
00252 edm::ParameterSet thePSet;
00253 LinearizationPointFinder* theLinP;
00254 VertexUpdator<N> * theUpdator;
00255 VertexSmoother<N> * theSmoother;
00256 const AbstractLTSFactory<N> * theLTrackFactory;
00257
00258 VertexTrackFactory<N> theVTrackFactory;
00259
00260
00261 };
00262
00263 #endif