CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SequentialVertexFitter.h
Go to the documentation of this file.
1 #ifndef SequentialVertexFitter_H
2 #define SequentialVertexFitter_H
3 
11 #include <cmath>
12 // #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
13 // #include "Vertex/VertexPrimitives/interface/VertexSeedFactory.h"
14 
33 template <unsigned int N>
35 
36 public:
37 
41 
42 
50  const VertexUpdator<N> & updator, const VertexSmoother<N> & smoother,
51  const AbstractLTSFactory<N> & ltsf);
52 
58  const LinearizationPointFinder & linP,
59  const VertexUpdator<N> & updator, const VertexSmoother<N> & smoother,
60  const AbstractLTSFactory<N> & ltsf);
61 
67 
68 
69  virtual ~SequentialVertexFitter();
70 
71 
79 
80 
85  void setMaximumNumberOfIterations(int maxIterations)
86  {theMaxStep = maxIterations;}
87 
95  virtual CachingVertex<N> vertex(const std::vector<reco::TransientTrack> & tracks) const;
96 
106  virtual CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack> & tracks) const;
107 
111  virtual CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack> & tracks, const reco::BeamSpot & spot ) const;
112 
113 
116  virtual CachingVertex<N> vertex(const std::vector<reco::TransientTrack> & tracks,
117  const GlobalPoint& linPoint) const;
118 
123  virtual CachingVertex<N> vertex(const std::vector<reco::TransientTrack> & tracks,
124  const reco::BeamSpot& beamSpot) const;
125 
126 
132  virtual CachingVertex<N> vertex(const std::vector<reco::TransientTrack> & tracks,
133  const GlobalPoint& priorPos,
134  const GlobalError& priorError) const;
135 
140  virtual CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack> & tracks,
141  const GlobalPoint& priorPos,
142  const GlobalError& priorError) const;
143 
144 
145 
156 // virtual CachingVertex<N> vertex(const RefCountedVertexSeed seed) const;
157 
162  {return theLinP;}
163 
165  {return theUpdator;}
166 
168  {return theSmoother;}
169 
170  const float maxShift() const
171  {return theMaxShift;}
172 
173  const int maxStep() const
174  {return theMaxStep;}
175 
177  {return thePSet;}
178 
180  return new SequentialVertexFitter(* this);
181  }
182 
184  { return theLTrackFactory;}
185 
186 protected:
187 
193 
194 private:
195 
205  CachingVertex<N> fit(const std::vector<RefCountedVertexTrack> & tracks,
206  const VertexState priorVertex, bool withPrior) const;
207 
215  std::vector<RefCountedVertexTrack> linearizeTracks(const std::vector<reco::TransientTrack> & tracks,
216  const VertexState state) const;
217 
228  std::vector<RefCountedVertexTrack> reLinearizeTracks(
229  const std::vector<RefCountedVertexTrack> & tracks,
230  const VertexState state) const;
231 
232 
237  void readParameters();
238  void setDefaultParameters();
239 
243  inline bool hasNan(const GlobalPoint& point) const {
244  using namespace std;
245  return (isnan(point.x())|| isnan(point.y()) || isnan(point.z()));
246  }
247 
248 
249  float theMaxShift;
251 
257 // LinearizedTrackStateFactoryAbstractLTSFactory theLTrackFactory;
259 
260  // VertexSeedFactory theVSeedFactory;
261 };
262 
263 #endif
std::vector< RefCountedVertexTrack > reLinearizeTracks(const std::vector< RefCountedVertexTrack > &tracks, const VertexState state) const
LinearizationPointFinder * theLinP
const VertexUpdator< N > * vertexUpdator() const
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
VertexUpdator< N > * theUpdator
list original
Definition: definitions.py:60
VertexTrackFactory< N > theVTrackFactory
T y() const
Definition: PV3DBase.h:62
const AbstractLTSFactory< N > * theLTrackFactory
const AbstractLTSFactory< N > * linearizedTrackStateFactory() const
virtual CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const
ReferenceCountingPointer< RefittedTrackState< N > > RefCountedRefittedTrackState
ReferenceCountingPointer< VertexTrack< N > > RefCountedVertexTrack
const VertexSmoother< N > * vertexSmoother() const
bool isnan(float x)
Definition: math.h:13
void setMaximumNumberOfIterations(int maxIterations)
T z() const
Definition: PV3DBase.h:63
SequentialVertexFitter * clone() const
const edm::ParameterSet parameterSet() const
tuple tracks
Definition: testEve_cfg.py:39
char state
Definition: procUtils.cc:75
bool hasNan(const GlobalPoint &point) const
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
T x() const
Definition: PV3DBase.h:61
ReferenceCountingPointer< LinearizedTrackState< N > > RefCountedLinearizedTrackState
const float maxShift() const
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
VertexSmoother< N > * theSmoother
const LinearizationPointFinder * linearizationPointFinder() const
void setMaximumDistance(float maxShift)