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
00012
00013
00032 template <unsigned int N>
00033 class SequentialVertexFitter : public VertexFitter<N> {
00034
00035 public:
00036
00037 typedef ReferenceCountingPointer<RefittedTrackState<N> > RefCountedRefittedTrackState;
00038 typedef ReferenceCountingPointer<VertexTrack<N> > RefCountedVertexTrack;
00039 typedef ReferenceCountingPointer<LinearizedTrackState<N> > RefCountedLinearizedTrackState;
00040
00041
00048 SequentialVertexFitter(const LinearizationPointFinder & linP,
00049 const VertexUpdator<N> & updator, const VertexSmoother<N> & smoother,
00050 const AbstractLTSFactory<N> & ltsf);
00051
00056 SequentialVertexFitter(const edm::ParameterSet& pSet,
00057 const LinearizationPointFinder & linP,
00058 const VertexUpdator<N> & updator, const VertexSmoother<N> & smoother,
00059 const AbstractLTSFactory<N> & ltsf);
00060
00065 SequentialVertexFitter(const SequentialVertexFitter & original);
00066
00067
00068 virtual ~SequentialVertexFitter();
00069
00070
00077 void setMaximumDistance(float maxShift) {theMaxShift = maxShift;}
00078
00079
00084 void setMaximumNumberOfIterations(int maxIterations)
00085 {theMaxStep = maxIterations;}
00086
00094 virtual CachingVertex<N> vertex(const std::vector<reco::TransientTrack> & tracks) const;
00095
00105 virtual CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack> & tracks) const;
00106
00110 virtual CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack> & tracks, const reco::BeamSpot & spot ) const;
00111
00112
00115 virtual CachingVertex<N> vertex(const std::vector<reco::TransientTrack> & tracks,
00116 const GlobalPoint& linPoint) const;
00117
00122 virtual CachingVertex<N> vertex(const vector<reco::TransientTrack> & tracks,
00123 const reco::BeamSpot& beamSpot) const;
00124
00125
00131 virtual CachingVertex<N> vertex(const std::vector<reco::TransientTrack> & tracks,
00132 const GlobalPoint& priorPos,
00133 const GlobalError& priorError) const;
00134
00139 virtual CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack> & tracks,
00140 const GlobalPoint& priorPos,
00141 const GlobalError& priorError) const;
00142
00143
00144
00155
00156
00160 const LinearizationPointFinder * linearizationPointFinder() const
00161 {return theLinP;}
00162
00163 const VertexUpdator<N> * vertexUpdator() const
00164 {return theUpdator;}
00165
00166 const VertexSmoother<N> * vertexSmoother() const
00167 {return theSmoother;}
00168
00169 const float maxShift() const
00170 {return theMaxShift;}
00171
00172 const int maxStep() const
00173 {return theMaxStep;}
00174
00175 const edm::ParameterSet parameterSet() const
00176 {return thePSet;}
00177
00178 SequentialVertexFitter * clone() const {
00179 return new SequentialVertexFitter(* this);
00180 }
00181
00182 const AbstractLTSFactory<N> * linearizedTrackStateFactory() const
00183 { return theLTrackFactory;}
00184
00185 protected:
00186
00191 SequentialVertexFitter() {}
00192
00193 private:
00194
00204 CachingVertex<N> fit(const std::vector<RefCountedVertexTrack> & tracks,
00205 const VertexState priorVertex, bool withPrior) const;
00206
00214 std::vector<RefCountedVertexTrack> linearizeTracks(const std::vector<reco::TransientTrack> & tracks,
00215 const VertexState state) const;
00216
00227 std::vector<RefCountedVertexTrack> reLinearizeTracks(
00228 const std::vector<RefCountedVertexTrack> & tracks,
00229 const VertexState state) const;
00230
00231
00236 void readParameters();
00237 void setDefaultParameters();
00238
00242 inline bool hasNan(const GlobalPoint& point) const {
00243 return (isnan(point.x())|| isnan(point.y()) || isnan(point.z()));
00244 }
00245
00246
00247 float theMaxShift;
00248 int theMaxStep;
00249
00250 edm::ParameterSet thePSet;
00251 LinearizationPointFinder* theLinP;
00252 VertexUpdator<N> * theUpdator;
00253 VertexSmoother<N> * theSmoother;
00254 const AbstractLTSFactory<N> * theLTrackFactory;
00255
00256 VertexTrackFactory<N> theVTrackFactory;
00257
00258
00259 };
00260
00261 #endif