CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Types | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
SequentialVertexFitter< N > Class Template Reference

#include <SequentialVertexFitter.h>

Inheritance diagram for SequentialVertexFitter< N >:
VertexFitter< N >

Public Types

typedef
ReferenceCountingPointer
< LinearizedTrackState< N > > 
RefCountedLinearizedTrackState
 
typedef
ReferenceCountingPointer
< RefittedTrackState< N > > 
RefCountedRefittedTrackState
 
typedef
ReferenceCountingPointer
< VertexTrack< N > > 
RefCountedVertexTrack
 

Public Member Functions

SequentialVertexFitterclone () const override
 
const LinearizationPointFinderlinearizationPointFinder () const
 
const AbstractLTSFactory< N > * linearizedTrackStateFactory () const
 
const float maxShift () const
 
const int maxStep () const
 
const edm::ParameterSet parameterSet () const
 
 SequentialVertexFitter (const LinearizationPointFinder &linP, const VertexUpdator< N > &updator, const VertexSmoother< N > &smoother, const AbstractLTSFactory< N > &ltsf)
 
 SequentialVertexFitter (const edm::ParameterSet &pSet, const LinearizationPointFinder &linP, const VertexUpdator< N > &updator, const VertexSmoother< N > &smoother, const AbstractLTSFactory< N > &ltsf)
 
 SequentialVertexFitter (const SequentialVertexFitter &original)
 
void setMaximumDistance (float maxShift)
 
void setMaximumNumberOfIterations (int maxIterations)
 
CachingVertex< Nvertex (const std::vector< reco::TransientTrack > &tracks) const override
 
CachingVertex< Nvertex (const std::vector< RefCountedVertexTrack > &tracks) const override
 
CachingVertex< Nvertex (const std::vector< RefCountedVertexTrack > &tracks, const reco::BeamSpot &spot) const override
 
CachingVertex< Nvertex (const std::vector< reco::TransientTrack > &tracks, const GlobalPoint &linPoint) const override
 
CachingVertex< Nvertex (const std::vector< reco::TransientTrack > &tracks, const reco::BeamSpot &beamSpot) const override
 
CachingVertex< Nvertex (const std::vector< reco::TransientTrack > &tracks, const GlobalPoint &priorPos, const GlobalError &priorError) const override
 
CachingVertex< Nvertex (const std::vector< RefCountedVertexTrack > &tracks, const GlobalPoint &priorPos, const GlobalError &priorError) const override
 
const VertexSmoother< N > * vertexSmoother () const
 
const VertexUpdator< N > * vertexUpdator () const
 
 ~SequentialVertexFitter () override
 
- Public Member Functions inherited from VertexFitter< N >
virtual CachingVertex< Nvertex (const std::vector< typename CachingVertex< N >::RefCountedVertexTrack > &tracks) const =0
 
virtual CachingVertex< Nvertex (const std::vector< typename CachingVertex< N >::RefCountedVertexTrack > &tracks, const reco::BeamSpot &spot) const =0
 
virtual CachingVertex< Nvertex (const std::vector< typename CachingVertex< N >::RefCountedVertexTrack > &tracks, const GlobalPoint &priorPos, const GlobalError &priorError) const =0
 
 VertexFitter ()
 
virtual ~VertexFitter ()
 

Protected Member Functions

 SequentialVertexFitter ()
 

Private Member Functions

CachingVertex< Nfit (const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
 
bool hasNan (const GlobalPoint &point) const
 
std::vector
< RefCountedVertexTrack
linearizeTracks (const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
 
void readParameters ()
 
std::vector
< RefCountedVertexTrack
reLinearizeTracks (const std::vector< RefCountedVertexTrack > &tracks, const VertexState state) const
 
void setDefaultParameters ()
 

Private Attributes

LinearizationPointFindertheLinP
 
const AbstractLTSFactory< N > * theLTrackFactory
 
float theMaxShift
 
int theMaxStep
 
edm::ParameterSet thePSet
 
VertexSmoother< N > * theSmoother
 
VertexUpdator< N > * theUpdator
 
VertexTrackFactory< NtheVTrackFactory
 

Detailed Description

template<unsigned int N>
class SequentialVertexFitter< N >

Sequential vertex fitter, where the vertex is updated one track at the time. The fitter will iterate over the set of tracks until the transverse distance between vertices computed in the previous and the current iterations is less than the specified convergence criteria, or until the maximum number of iterations is reached. The transverse distance determines the linearization error. The default convergence criterion is 1 mm. The default maximum number of steps is 10. These parameters can be configured in .orcarc ( SequentialVertexFitter:maximumDistance and SequentialVertexFitter:maximumNumberOfIterations). After the vertex fit, the tracks can be refit with the additional constraint of the vertex position.

Definition at line 34 of file SequentialVertexFitter.h.

Member Typedef Documentation

Definition at line 38 of file SequentialVertexFitter.h.

Definition at line 36 of file SequentialVertexFitter.h.

Definition at line 37 of file SequentialVertexFitter.h.

Constructor & Destructor Documentation

template<unsigned int N>
SequentialVertexFitter< N >::SequentialVertexFitter ( const LinearizationPointFinder linP,
const VertexUpdator< N > &  updator,
const VertexSmoother< N > &  smoother,
const AbstractLTSFactory< N > &  ltsf 
)

Reimplemented constructors to use any kind of linearisation point finder, vertex updator and smoother. If no smoother is to be used, do not specify an instance for it.

Definition at line 21 of file SequentialVertexFitter.cc.

References SequentialVertexFitter< N >::setDefaultParameters().

25  : theLinP(linP.clone()),
26  theUpdator(updator.clone()),
27  theSmoother(smoother.clone()),
28  theLTrackFactory(ltsf.clone()) {
30 }
LinearizationPointFinder * theLinP
VertexUpdator< N > * theUpdator
virtual VertexSmoother * clone() const =0
const AbstractLTSFactory< N > * theLTrackFactory
virtual const AbstractLTSFactory * clone() const =0
virtual LinearizationPointFinder * clone() const =0
virtual VertexUpdator * clone() const =0
VertexSmoother< N > * theSmoother
template<unsigned int N>
SequentialVertexFitter< N >::SequentialVertexFitter ( const edm::ParameterSet pSet,
const LinearizationPointFinder linP,
const VertexUpdator< N > &  updator,
const VertexSmoother< N > &  smoother,
const AbstractLTSFactory< N > &  ltsf 
)

Same as above, using a ParameterSet to set the convergence criteria

Definition at line 33 of file SequentialVertexFitter.cc.

References SequentialVertexFitter< N >::readParameters().

38  : thePSet(pSet),
39  theLinP(linP.clone()),
40  theUpdator(updator.clone()),
41  theSmoother(smoother.clone()),
42  theLTrackFactory(ltsf.clone()) {
44 }
LinearizationPointFinder * theLinP
VertexUpdator< N > * theUpdator
virtual VertexSmoother * clone() const =0
const AbstractLTSFactory< N > * theLTrackFactory
virtual const AbstractLTSFactory * clone() const =0
virtual LinearizationPointFinder * clone() const =0
virtual VertexUpdator * clone() const =0
VertexSmoother< N > * theSmoother
template<unsigned int N>
SequentialVertexFitter< N >::SequentialVertexFitter ( const SequentialVertexFitter< N > &  original)

Copy constructor

Definition at line 47 of file SequentialVertexFitter.cc.

References LinearizationPointFinder::clone(), SequentialVertexFitter< N >::linearizationPointFinder(), SequentialVertexFitter< N >::linearizedTrackStateFactory(), SequentialVertexFitter< N >::maxShift(), SequentialVertexFitter< N >::maxStep(), SequentialVertexFitter< N >::parameterSet(), SequentialVertexFitter< N >::vertexSmoother(), and SequentialVertexFitter< N >::vertexUpdator().

47  {
48  thePSet = original.parameterSet();
49  theLinP = original.linearizationPointFinder()->clone();
50  theUpdator = original.vertexUpdator()->clone();
51  theSmoother = original.vertexSmoother()->clone();
52  theMaxShift = original.maxShift();
53  theMaxStep = original.maxStep();
54  theLTrackFactory = original.linearizedTrackStateFactory()->clone();
55 }
LinearizationPointFinder * theLinP
const VertexUpdator< N > * vertexUpdator() const
VertexUpdator< N > * theUpdator
const AbstractLTSFactory< N > * theLTrackFactory
const AbstractLTSFactory< N > * linearizedTrackStateFactory() const
const VertexSmoother< N > * vertexSmoother() const
const edm::ParameterSet parameterSet() const
virtual LinearizationPointFinder * clone() const =0
const float maxShift() const
VertexSmoother< N > * theSmoother
const LinearizationPointFinder * linearizationPointFinder() const
template<unsigned int N>
SequentialVertexFitter< N >::~SequentialVertexFitter ( )
override

Definition at line 58 of file SequentialVertexFitter.cc.

58  {
59  delete theLinP;
60  delete theUpdator;
61  delete theSmoother;
62  delete theLTrackFactory;
63 }
LinearizationPointFinder * theLinP
VertexUpdator< N > * theUpdator
const AbstractLTSFactory< N > * theLTrackFactory
VertexSmoother< N > * theSmoother
template<unsigned int N>
SequentialVertexFitter< N >::SequentialVertexFitter ( )
inlineprotected

Default constructor. Is here, as we do not want anybody to use it.

Definition at line 172 of file SequentialVertexFitter.h.

Referenced by SequentialVertexFitter< 5 >::clone().

172 {}

Member Function Documentation

template<unsigned int N>
SequentialVertexFitter* SequentialVertexFitter< N >::clone ( ) const
inlineoverridevirtual

Fit vertex out of a VertexSeed

Implements VertexFitter< N >.

Definition at line 163 of file SequentialVertexFitter.h.

Referenced by GsfVertexFitter::GsfVertexFitter().

163 { return new SequentialVertexFitter(*this); }
template<unsigned int N>
CachingVertex< N > SequentialVertexFitter< N >::fit ( const std::vector< RefCountedVertexTrack > &  tracks,
const VertexState  priorVertex,
bool  withPrior 
) const
private

The methode where the vrte fit is actually done. The seed is used as the prior estimate in the vertex fit (in case its error is large, it will have little influence on the fit. The tracks will be relinearized in case further loops are needed. tracks The tracks to use in the fit. priorVertex The prior estimate of the vertex

Returns
The fitted vertex

Definition at line 221 of file SequentialVertexFitter.cc.

References relativeConstraints::error, VertexState::error(), mps_fire::i, gpuClustering::id, CachingVertex< N >::isValid(), LogDebug, VertexState::position(), CachingVertex< N >::position(), tracks, transverse(), and CachingVertex< N >::vertexState().

Referenced by trackingPlots.Iteration::modules().

223  {
224  std::vector<RefCountedVertexTrack> initialTracks;
225  GlobalPoint priorVertexPosition = priorVertex.position();
226  GlobalError priorVertexError = priorVertex.error();
227 
228  CachingVertex<N> returnVertex(priorVertexPosition, priorVertexError, initialTracks, 0);
229  if (withPrior) {
230  returnVertex = CachingVertex<N>(
231  priorVertexPosition, priorVertexError, priorVertexPosition, priorVertexError, initialTracks, 0);
232  }
233  CachingVertex<N> initialVertex = returnVertex;
234  std::vector<RefCountedVertexTrack> globalVTracks = tracks;
235 
236  // main loop through all the VTracks
237  bool validVertex = true;
238  int step = 0;
239  GlobalPoint newPosition = priorVertexPosition;
240  GlobalPoint previousPosition;
241  do {
242  CachingVertex<N> fVertex = initialVertex;
243  // make new linearized and vertex tracks for the next iteration
244  if (step != 0)
245  globalVTracks = reLinearizeTracks(tracks, returnVertex.vertexState());
246 
247  // update sequentially the vertex estimate
248  for (typename std::vector<RefCountedVertexTrack>::const_iterator i = globalVTracks.begin();
249  i != globalVTracks.end();
250  i++) {
251  fVertex = theUpdator->add(fVertex, *i);
252  if (!fVertex.isValid())
253  break;
254  }
255 
256  validVertex = fVertex.isValid();
257  // check tracker bounds and NaN in position
258  if (validVertex && hasNan(fVertex.position())) {
259  LogDebug("RecoVertex/SequentialVertexFitter") << "Fitted position is NaN.\n";
260  validVertex = false;
261  }
262 
263  if (validVertex && !insideTrackerBounds(fVertex.position())) {
264  LogDebug("RecoVertex/SequentialVertexFitter") << "Fitted position is out of tracker bounds.\n";
265  validVertex = false;
266  }
267 
268  if (!validVertex) {
269  // reset initial vertex position to (0,0,0) and force new iteration
270  // if number of steps not exceeded
271  ROOT::Math::SMatrixIdentity id;
272  AlgebraicSymMatrix33 we(id);
273  GlobalError error(we * 10000);
274  fVertex = CachingVertex<N>(GlobalPoint(0, 0, 0), error, initialTracks, 0);
275  }
276 
277  previousPosition = newPosition;
278  newPosition = fVertex.position();
279 
280  returnVertex = fVertex;
281  globalVTracks.clear();
282  step++;
283  } while ((step != theMaxStep) && (((previousPosition - newPosition).transverse() > theMaxShift) || (!validVertex)));
284 
285  if (!validVertex) {
286  LogDebug("RecoVertex/SequentialVertexFitter")
287  << "Fitted position is invalid (out of tracker bounds or has NaN). Returned vertex is invalid\n";
288  return CachingVertex<N>(); // return invalid vertex
289  }
290 
291  if (step >= theMaxStep) {
292  LogDebug("RecoVertex/SequentialVertexFitter")
293  << "The maximum number of steps has been exceeded. Returned vertex is invalid\n";
294  return CachingVertex<N>(); // return invalid vertex
295  }
296 
297  // smoothing
298  returnVertex = theSmoother->smooth(returnVertex);
299 
300  return returnVertex;
301 }
std::vector< RefCountedVertexTrack > reLinearizeTracks(const std::vector< RefCountedVertexTrack > &tracks, const VertexState state) const
uint16_t *__restrict__ id
VertexUpdator< N > * theUpdator
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
auto const & tracks
cannot be loose
GlobalPoint position() const
Definition: VertexState.h:62
GlobalPoint position() const
T transverse() const
Another name for perp()
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
bool isValid() const
GlobalError error() const
Definition: VertexState.h:64
bool hasNan(const GlobalPoint &point) const
step
Definition: StallMonitor.cc:98
VertexSmoother< N > * theSmoother
#define LogDebug(id)
template<unsigned int N>
bool SequentialVertexFitter< N >::hasNan ( const GlobalPoint point) const
inlineprivate

Checks whether any of the three coordinates is a Nan

Definition at line 221 of file SequentialVertexFitter.h.

221  {
222  using namespace std;
223  return (edm::isNotFinite(point.x()) || edm::isNotFinite(point.y()) || edm::isNotFinite(point.z()));
224  }
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
T y() const
Definition: PV3DBase.h:60
T z() const
Definition: PV3DBase.h:61
T x() const
Definition: PV3DBase.h:59
template<unsigned int N>
const LinearizationPointFinder* SequentialVertexFitter< N >::linearizationPointFinder ( ) const
inline

Method returning the fitted vertex, from a VertexSeed. For the first loop, the position of the VertexSeed will be used as linearization point. If subsequent loops are needed, the new VertexTracks will be created with the last estimate of the vertex as linearization point. In case a non-sero error is given, the position and error of the VertexSeed will be used as prior estimate in the vertex fit.

Parameters
seedThe VertexSeed to fit.
Returns
The fitted vertex Access methods

Definition at line 151 of file SequentialVertexFitter.h.

Referenced by SequentialVertexFitter< N >::SequentialVertexFitter().

151 { return theLinP; }
LinearizationPointFinder * theLinP
template<unsigned int N>
const AbstractLTSFactory<N>* SequentialVertexFitter< N >::linearizedTrackStateFactory ( ) const
inline

Definition at line 165 of file SequentialVertexFitter.h.

Referenced by SequentialVertexFitter< N >::SequentialVertexFitter().

165 { return theLTrackFactory; }
const AbstractLTSFactory< N > * theLTrackFactory
template<unsigned int N>
vector< typename SequentialVertexFitter< N >::RefCountedVertexTrack > SequentialVertexFitter< N >::linearizeTracks ( const std::vector< reco::TransientTrack > &  tracks,
const VertexState  state 
) const
private

Construct a container of VertexTrack from a set of RecTracks.

Parameters
tracksThe container of RecTracks.
seedThe seed to use for the VertexTracks. This position will also be used as the new linearization point.
Returns
The container of VertexTracks which are to be used in the next fit.

Definition at line 184 of file SequentialVertexFitter.cc.

References mps_fire::i, and VertexState::position().

185  {
186  GlobalPoint linP = state.position();
187  std::vector<RefCountedVertexTrack> finalTracks;
188  finalTracks.reserve(tracks.size());
189  for (vector<reco::TransientTrack>::const_iterator i = tracks.begin(); i != tracks.end(); i++) {
190  RefCountedLinearizedTrackState lTrData = theLTrackFactory->linearizedTrackState(linP, *i);
191  RefCountedVertexTrack vTrData = theVTrackFactory.vertexTrack(lTrData, state);
192  finalTracks.push_back(vTrData);
193  }
194  return finalTracks;
195 }
VertexTrackFactory< N > theVTrackFactory
const AbstractLTSFactory< N > * theLTrackFactory
auto const & tracks
cannot be loose
GlobalPoint position() const
Definition: VertexState.h:62
ReferenceCountingPointer< VertexTrack< N > > RefCountedVertexTrack
ReferenceCountingPointer< LinearizedTrackState< N > > RefCountedLinearizedTrackState
template<unsigned int N>
const float SequentialVertexFitter< N >::maxShift ( ) const
inline
template<unsigned int N>
const int SequentialVertexFitter< N >::maxStep ( ) const
inline
template<unsigned int N>
const edm::ParameterSet SequentialVertexFitter< N >::parameterSet ( ) const
inline

Definition at line 161 of file SequentialVertexFitter.h.

Referenced by SequentialVertexFitter< N >::SequentialVertexFitter().

161 { return thePSet; }
template<unsigned int N>
void SequentialVertexFitter< N >::readParameters ( )
private

Reads the configurable parameters.

Definition at line 66 of file SequentialVertexFitter.cc.

Referenced by SequentialVertexFitter< N >::SequentialVertexFitter().

66  {
67  theMaxShift = thePSet.getParameter<double>("maxDistance"); //0.01
68  theMaxStep = thePSet.getParameter<int>("maxNbrOfIterations"); //10
69 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
template<unsigned int N>
vector< typename SequentialVertexFitter< N >::RefCountedVertexTrack > SequentialVertexFitter< N >::reLinearizeTracks ( const std::vector< RefCountedVertexTrack > &  tracks,
const VertexState  state 
) const
private

Construct new a container of VertexTrack with a new linearization point and vertex seed, from an existing set of VertexTrack, from which only the recTracks will be used.

Parameters
tracksThe original container of VertexTracks, from which the RecTracks will be extracted.
seedThe seed to use for the VertexTracks. This position will also be used as the new linearization point.
Returns
The container of VertexTracks which are to be used in the next fit.

Definition at line 202 of file SequentialVertexFitter.cc.

References mps_fire::i, VertexState::position(), and VertexState::weight().

203  {
204  GlobalPoint linP = state.position();
205  std::vector<RefCountedVertexTrack> finalTracks;
206  finalTracks.reserve(tracks.size());
207  for (typename std::vector<RefCountedVertexTrack>::const_iterator i = tracks.begin(); i != tracks.end(); i++) {
208  RefCountedLinearizedTrackState lTrData = (**i).linearizedTrack()->stateWithNewLinearizationPoint(linP);
209  // RefCountedLinearizedTrackState lTrData =
210  // theLTrackFactory->linearizedTrackState(linP,
211  // (**i).linearizedTrack()->track());
212  RefCountedVertexTrack vTrData = theVTrackFactory.vertexTrack(lTrData, state, (**i).weight());
213  finalTracks.push_back(vTrData);
214  }
215  return finalTracks;
216 }
VertexTrackFactory< N > theVTrackFactory
auto const & tracks
cannot be loose
GlobalPoint position() const
Definition: VertexState.h:62
ReferenceCountingPointer< VertexTrack< N > > RefCountedVertexTrack
GlobalWeight weight() const
Definition: VertexState.h:69
ReferenceCountingPointer< LinearizedTrackState< N > > RefCountedLinearizedTrackState
template<unsigned int N>
void SequentialVertexFitter< N >::setDefaultParameters ( )
private

Definition at line 72 of file SequentialVertexFitter.cc.

Referenced by SequentialVertexFitter< N >::SequentialVertexFitter().

72  {
73  thePSet.addParameter<double>("maxDistance", 0.01);
74  thePSet.addParameter<int>("maxNbrOfIterations", 10); //10
76 }
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:135
template<unsigned int N>
void SequentialVertexFitter< N >::setMaximumDistance ( float  maxShift)
inline

Method to set the convergence criterion (the maximum distance between the vertex computed in the previous and the current iterations to consider the fit to have converged)

Definition at line 75 of file SequentialVertexFitter.h.

template<unsigned int N>
void SequentialVertexFitter< N >::setMaximumNumberOfIterations ( int  maxIterations)
inline

Method to set the maximum number of iterations to perform

Definition at line 81 of file SequentialVertexFitter.h.

81 { theMaxStep = maxIterations; }
template<unsigned int N>
CachingVertex< N > SequentialVertexFitter< N >::vertex ( const std::vector< reco::TransientTrack > &  tracks) const
overridevirtual

Method returning the fitted vertex, from a container of reco::TransientTracks. The linearization point will be searched with the given LP finder. No prior vertex position will be used in the vertex fit.

Parameters
tracksThe container of RecTracks to fit.
Returns
The fitted vertex

Implements VertexFitter< N >.

Definition at line 79 of file SequentialVertexFitter.cc.

References relativeConstraints::error, gpuVertexFinder::fit, and gpuClustering::id.

Referenced by Tau.Tau::dxy(), GsfVertexFitter::vertex(), and KalmanVertexFitter::vertex().

79  {
80  // Linearization Point
82  if (!insideTrackerBounds(linP))
83  linP = GlobalPoint(0, 0, 0);
84 
85  // Initial vertex state, with a very large error matrix
86  ROOT::Math::SMatrixIdentity id;
87  AlgebraicSymMatrix33 we(id);
88  GlobalError error(we * 10000);
89  VertexState state(linP, error);
90  std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
91  return fit(vtContainer, state, false);
92 }
uint16_t *__restrict__ id
LinearizationPointFinder * theLinP
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
auto const & tracks
cannot be loose
virtual GlobalPoint getLinearizationPoint(const std::vector< reco::TransientTrack > &) const =0
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
template<unsigned int N>
CachingVertex< N > SequentialVertexFitter< N >::vertex ( const std::vector< RefCountedVertexTrack > &  tracks) const
override

Method returning the fitted vertex, from a container of VertexTracks. For the first loop, the LinearizedTrackState contained in the VertexTracks will be used. If subsequent loops are needed, the new VertexTracks will be created with the last estimate of the vertex as linearization point. No prior vertex position will be used in the vertex fit.

Parameters
tracksThe container of VertexTracks to fit.
Returns
The fitted vertex

Definition at line 102 of file SequentialVertexFitter.cc.

References relativeConstraints::error, gpuVertexFinder::fit, and gpuClustering::id.

Referenced by Tau.Tau::dxy().

102  {
103  // Initial vertex state, with a very small weight matrix
104  GlobalPoint linP = tracks[0]->linearizedTrack()->linearizationPoint();
105  ROOT::Math::SMatrixIdentity id;
106  AlgebraicSymMatrix33 we(id);
107  GlobalError error(we * 10000);
108  VertexState state(linP, error);
109  return fit(tracks, state, false);
110 }
uint16_t *__restrict__ id
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
auto const & tracks
cannot be loose
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
template<unsigned int N>
CachingVertex< N > SequentialVertexFitter< N >::vertex ( const std::vector< RefCountedVertexTrack > &  tracks,
const reco::BeamSpot spot 
) const
override

Same as above, only now also with BeamSpot!

Definition at line 95 of file SequentialVertexFitter.cc.

References gpuVertexFinder::fit.

Referenced by Tau.Tau::dxy().

96  {
97  VertexState state(spot);
98  return fit(tracks, state, true);
99 }
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
auto const & tracks
cannot be loose
template<unsigned int N>
CachingVertex< N > SequentialVertexFitter< N >::vertex ( const std::vector< reco::TransientTrack > &  tracks,
const GlobalPoint linPoint 
) const
overridevirtual

Fit vertex out of a set of RecTracks. Uses the specified linearization point.

Implements VertexFitter< N >.

Definition at line 116 of file SequentialVertexFitter.cc.

References relativeConstraints::error, gpuVertexFinder::fit, and gpuClustering::id.

Referenced by Tau.Tau::dxy().

117  {
118  // Initial vertex state, with a very large error matrix
119  ROOT::Math::SMatrixIdentity id;
120  AlgebraicSymMatrix33 we(id);
121  GlobalError error(we * 10000);
122  VertexState state(linPoint, error);
123  std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
124  return fit(vtContainer, state, false);
125 }
uint16_t *__restrict__ id
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
auto const & tracks
cannot be loose
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
template<unsigned int N>
CachingVertex< N > SequentialVertexFitter< N >::vertex ( const std::vector< reco::TransientTrack > &  tracks,
const reco::BeamSpot beamSpot 
) const
overridevirtual

Fit vertex out of a set of TransientTracks. The specified BeamSpot will be used as priot, but NOT for the linearization. The specified LinearizationPointFinder will be used to find the linearization point.

Implements VertexFitter< N >.

Definition at line 132 of file SequentialVertexFitter.cc.

References relativeConstraints::error, gpuVertexFinder::fit, and gpuClustering::id.

Referenced by Tau.Tau::dxy().

133  {
134  VertexState beamSpotState(beamSpot);
135  std::vector<RefCountedVertexTrack> vtContainer;
136 
137  if (tracks.size() > 1) {
138  // Linearization Point search if there are more than 1 track
140  if (!insideTrackerBounds(linP))
141  linP = GlobalPoint(0, 0, 0);
142  ROOT::Math::SMatrixIdentity id;
143  AlgebraicSymMatrix33 we(id);
144  GlobalError error(we * 10000);
145  VertexState lpState(linP, error);
146  vtContainer = linearizeTracks(tracks, lpState);
147  } else {
148  // otherwise take the beamspot position.
149  vtContainer = linearizeTracks(tracks, beamSpotState);
150  }
151 
152  return fit(vtContainer, beamSpotState, true);
153 }
uint16_t *__restrict__ id
LinearizationPointFinder * theLinP
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
auto const & tracks
cannot be loose
virtual GlobalPoint getLinearizationPoint(const std::vector< reco::TransientTrack > &) const =0
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
template<unsigned int N>
CachingVertex< N > SequentialVertexFitter< N >::vertex ( const std::vector< reco::TransientTrack > &  tracks,
const GlobalPoint priorPos,
const GlobalError priorError 
) const
overridevirtual

Fit vertex out of a set of RecTracks. Uses the position as both the linearization point AND as prior estimate of the vertex position. The error is used for the weight of the prior estimate.

Implements VertexFitter< N >.

Definition at line 161 of file SequentialVertexFitter.cc.

References gpuVertexFinder::fit.

Referenced by Tau.Tau::dxy().

163  {
164  VertexState state(priorPos, priorError);
165  std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
166  return fit(vtContainer, state, true);
167 }
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
auto const & tracks
cannot be loose
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
template<unsigned int N>
CachingVertex< N > SequentialVertexFitter< N >::vertex ( const std::vector< RefCountedVertexTrack > &  tracks,
const GlobalPoint priorPos,
const GlobalError priorError 
) const
override

Fit vertex out of a set of VertexTracks Uses the position and error for the prior estimate of the vertex. This position is not used to relinearize the tracks.

Definition at line 174 of file SequentialVertexFitter.cc.

References gpuVertexFinder::fit.

Referenced by Tau.Tau::dxy().

176  {
177  VertexState state(priorPos, priorError);
178  return fit(tracks, state, true);
179 }
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
auto const & tracks
cannot be loose
template<unsigned int N>
const VertexSmoother<N>* SequentialVertexFitter< N >::vertexSmoother ( ) const
inline

Definition at line 155 of file SequentialVertexFitter.h.

Referenced by SequentialVertexFitter< N >::SequentialVertexFitter().

155 { return theSmoother; }
VertexSmoother< N > * theSmoother
template<unsigned int N>
const VertexUpdator<N>* SequentialVertexFitter< N >::vertexUpdator ( ) const
inline

Definition at line 153 of file SequentialVertexFitter.h.

Referenced by SequentialVertexFitter< N >::SequentialVertexFitter().

153 { return theUpdator; }
VertexUpdator< N > * theUpdator

Member Data Documentation

template<unsigned int N>
LinearizationPointFinder* SequentialVertexFitter< N >::theLinP
private
template<unsigned int N>
const AbstractLTSFactory<N>* SequentialVertexFitter< N >::theLTrackFactory
private
template<unsigned int N>
float SequentialVertexFitter< N >::theMaxShift
private
template<unsigned int N>
int SequentialVertexFitter< N >::theMaxStep
private
template<unsigned int N>
edm::ParameterSet SequentialVertexFitter< N >::thePSet
private

Definition at line 229 of file SequentialVertexFitter.h.

Referenced by SequentialVertexFitter< 5 >::parameterSet().

template<unsigned int N>
VertexSmoother<N>* SequentialVertexFitter< N >::theSmoother
private
template<unsigned int N>
VertexUpdator<N>* SequentialVertexFitter< N >::theUpdator
private
template<unsigned int N>
VertexTrackFactory<N> SequentialVertexFitter< N >::theVTrackFactory
private

Definition at line 235 of file SequentialVertexFitter.h.