CMS 3D CMS Logo

SequentialVertexFitter.h
Go to the documentation of this file.
1 #ifndef SequentialVertexFitter_H
2 #define SequentialVertexFitter_H
3 
12 #include <cmath>
13 // #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
14 // #include "Vertex/VertexPrimitives/interface/VertexSeedFactory.h"
15 
33 template <unsigned int N>
35 public:
39 
48  const VertexSmoother<N>& smoother,
49  const AbstractLTSFactory<N>& ltsf);
50 
56  const LinearizationPointFinder& linP,
58  const VertexSmoother<N>& smoother,
59  const AbstractLTSFactory<N>& ltsf);
60 
66 
67  ~SequentialVertexFitter() override;
68 
76 
82 
90  CachingVertex<N> vertex(const std::vector<reco::TransientTrack>& tracks) const override;
91 
101  CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack>& tracks) const override;
102 
106  CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack>& tracks, const reco::BeamSpot& spot) const override;
107 
110  CachingVertex<N> vertex(const std::vector<reco::TransientTrack>& tracks, const GlobalPoint& linPoint) const override;
111 
116  CachingVertex<N> vertex(const std::vector<reco::TransientTrack>& tracks,
117  const reco::BeamSpot& beamSpot) const override;
118 
124  CachingVertex<N> vertex(const std::vector<reco::TransientTrack>& tracks,
125  const GlobalPoint& priorPos,
126  const GlobalError& priorError) const override;
127 
132  CachingVertex<N> vertex(const std::vector<RefCountedVertexTrack>& tracks,
133  const GlobalPoint& priorPos,
134  const GlobalError& priorError) const override;
135 
146  // virtual CachingVertex<N> vertex(const RefCountedVertexSeed seed) const;
147 
152 
153  const VertexUpdator<N>* vertexUpdator() const { return theUpdator; }
154 
155  const VertexSmoother<N>* vertexSmoother() const { return theSmoother; }
156 
157  const float maxShift() const { return theMaxShift; }
158 
159  const int maxStep() const { return theMaxStep; }
160 
161  const edm::ParameterSet parameterSet() const { return thePSet; }
162 
163  SequentialVertexFitter* clone() const override { return new SequentialVertexFitter(*this); }
164 
166 
167 protected:
173 
174 private:
184  CachingVertex<N> fit(const std::vector<RefCountedVertexTrack>& tracks,
185  const VertexState priorVertex,
186  bool withPrior) const;
187 
195  std::vector<RefCountedVertexTrack> linearizeTracks(const std::vector<reco::TransientTrack>& tracks,
196  const VertexState state) const;
197 
208  std::vector<RefCountedVertexTrack> reLinearizeTracks(const std::vector<RefCountedVertexTrack>& tracks,
209  const VertexState state) const;
210 
215  void readParameters();
216  void setDefaultParameters();
217 
221  inline bool hasNan(const GlobalPoint& point) const {
222  using namespace std;
223  return (edm::isNotFinite(point.x()) || edm::isNotFinite(point.y()) || edm::isNotFinite(point.z()));
224  }
225 
226  float theMaxShift;
228 
234  // LinearizedTrackStateFactoryAbstractLTSFactory theLTrackFactory;
236 
237  // VertexSeedFactory theVSeedFactory;
238 };
239 
240 #endif
bool hasNan(const GlobalPoint &point) const
const AbstractLTSFactory< N > * linearizedTrackStateFactory() const
LinearizationPointFinder * theLinP
SequentialVertexFitter * clone() const override
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
VertexUpdator< N > * theUpdator
VertexTrackFactory< N > theVTrackFactory
const AbstractLTSFactory< N > * theLTrackFactory
std::vector< RefCountedVertexTrack > reLinearizeTracks(const std::vector< RefCountedVertexTrack > &tracks, const VertexState state) const
const VertexSmoother< N > * vertexSmoother() const
const edm::ParameterSet parameterSet() const
CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const override
const LinearizationPointFinder * linearizationPointFinder() const
ReferenceCountingPointer< RefittedTrackState< N > > RefCountedRefittedTrackState
ReferenceCountingPointer< VertexTrack< N > > RefCountedVertexTrack
void setMaximumNumberOfIterations(int maxIterations)
auto const & tracks
cannot be loose
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
const float maxShift() const
const VertexUpdator< N > * vertexUpdator() const
ReferenceCountingPointer< LinearizedTrackState< N > > RefCountedLinearizedTrackState
*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
void setMaximumDistance(float maxShift)