CMS 3D CMS Logo

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 bool insideTrackerBounds (const GlobalPoint &point) const
 
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)
 
void setTrackerBounds (float radius, float halfLength)
 
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< RefCountedVertexTracklinearizeTracks (const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
 
void readParameters ()
 
std::vector< RefCountedVertexTrackreLinearizeTracks (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
 
float trackerBoundsHalfLength {273.5}
 
float trackerBoundsRadius {112.}
 

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

◆ RefCountedLinearizedTrackState

Definition at line 38 of file SequentialVertexFitter.h.

◆ RefCountedRefittedTrackState

Definition at line 36 of file SequentialVertexFitter.h.

◆ RefCountedVertexTrack

Definition at line 37 of file SequentialVertexFitter.h.

Constructor & Destructor Documentation

◆ SequentialVertexFitter() [1/4]

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 10 of file SequentialVertexFitter.cc.

14  : theLinP(linP.clone()),
15  theUpdator(updator.clone()),
16  theSmoother(smoother.clone()),
17  theLTrackFactory(ltsf.clone()) {
19 }
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
VertexSmoother< N > * theSmoother

◆ SequentialVertexFitter() [2/4]

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 22 of file SequentialVertexFitter.cc.

27  : thePSet(pSet),
28  theLinP(linP.clone()),
29  theUpdator(updator.clone()),
30  theSmoother(smoother.clone()),
31  theLTrackFactory(ltsf.clone()) {
33 }
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
VertexSmoother< N > * theSmoother

◆ SequentialVertexFitter() [3/4]

template<unsigned int N>
SequentialVertexFitter< N >::SequentialVertexFitter ( const SequentialVertexFitter< N > &  original)

Copy constructor

Definition at line 36 of file SequentialVertexFitter.cc.

36  {
37  thePSet = original.parameterSet();
38  theLinP = original.linearizationPointFinder()->clone();
39  theUpdator = original.vertexUpdator()->clone();
40  theSmoother = original.vertexSmoother()->clone();
41  theMaxShift = original.maxShift();
42  theMaxStep = original.maxStep();
43  theLTrackFactory = original.linearizedTrackStateFactory()->clone();
44 }
LinearizationPointFinder * theLinP
VertexUpdator< N > * theUpdator
const AbstractLTSFactory< N > * theLTrackFactory
VertexSmoother< N > * theSmoother

◆ ~SequentialVertexFitter()

template<unsigned int N>
SequentialVertexFitter< N >::~SequentialVertexFitter ( )
override

Definition at line 47 of file SequentialVertexFitter.cc.

47  {
48  delete theLinP;
49  delete theUpdator;
50  delete theSmoother;
51  delete theLTrackFactory;
52 }
LinearizationPointFinder * theLinP
VertexUpdator< N > * theUpdator
const AbstractLTSFactory< N > * theLTrackFactory
VertexSmoother< N > * theSmoother

◆ SequentialVertexFitter() [4/4]

template<unsigned int N>
SequentialVertexFitter< N >::SequentialVertexFitter ( )
inlineprotected

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

Definition at line 188 of file SequentialVertexFitter.h.

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

188 {}

Member Function Documentation

◆ clone()

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.

163 { return new SequentialVertexFitter(*this); }

◆ fit()

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 210 of file SequentialVertexFitter.cc.

Referenced by trackingPlots.Iteration::modules().

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

◆ hasNan()

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 237 of file SequentialVertexFitter.h.

237  {
238  using namespace std;
239  return (edm::isNotFinite(point.x()) || edm::isNotFinite(point.y()) || edm::isNotFinite(point.z()));
240  }
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
*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

◆ insideTrackerBounds()

template<unsigned int N>
const bool SequentialVertexFitter< N >::insideTrackerBounds ( const GlobalPoint point) const
inline

Method checking whether a point is within the "tracker" bounds. The default values are set to the CMS inner tracker and vertices outsides these bounds will be rejected. To reconstruct vertices within the full detector including the muon system, set the tracker bounds to larger values.

Definition at line 174 of file SequentialVertexFitter.h.

174  {
175  return ((point.transverse() < trackerBoundsRadius) && (abs(point.z()) < trackerBoundsHalfLength));
176  }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
*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

◆ linearizationPointFinder()

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.

151 { return theLinP; }
LinearizationPointFinder * theLinP

◆ linearizedTrackStateFactory()

template<unsigned int N>
const AbstractLTSFactory<N>* SequentialVertexFitter< N >::linearizedTrackStateFactory ( ) const
inline

Definition at line 165 of file SequentialVertexFitter.h.

165 { return theLTrackFactory; }
const AbstractLTSFactory< N > * theLTrackFactory

◆ linearizeTracks()

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 173 of file SequentialVertexFitter.cc.

174  {
175  GlobalPoint linP = state.position();
176  std::vector<RefCountedVertexTrack> finalTracks;
177  finalTracks.reserve(tracks.size());
178  for (vector<reco::TransientTrack>::const_iterator i = tracks.begin(); i != tracks.end(); i++) {
179  RefCountedLinearizedTrackState lTrData = theLTrackFactory->linearizedTrackState(linP, *i);
180  RefCountedVertexTrack vTrData = theVTrackFactory.vertexTrack(lTrData, state);
181  finalTracks.push_back(vTrData);
182  }
183  return finalTracks;
184 }
VertexTrackFactory< N > theVTrackFactory
const AbstractLTSFactory< N > * theLTrackFactory
ReferenceCountingPointer< VertexTrack< N > > RefCountedVertexTrack
ReferenceCountingPointer< LinearizedTrackState< N > > RefCountedLinearizedTrackState

◆ maxShift()

template<unsigned int N>
const float SequentialVertexFitter< N >::maxShift ( ) const
inline

◆ maxStep()

template<unsigned int N>
const int SequentialVertexFitter< N >::maxStep ( ) const
inline

Definition at line 159 of file SequentialVertexFitter.h.

◆ parameterSet()

template<unsigned int N>
const edm::ParameterSet SequentialVertexFitter< N >::parameterSet ( ) const
inline

Definition at line 161 of file SequentialVertexFitter.h.

161 { return thePSet; }

◆ readParameters()

template<unsigned int N>
void SequentialVertexFitter< N >::readParameters ( )
private

Reads the configurable parameters.

Definition at line 55 of file SequentialVertexFitter.cc.

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

55  {
56  theMaxShift = thePSet.getParameter<double>("maxDistance"); //0.01
57  theMaxStep = thePSet.getParameter<int>("maxNbrOfIterations"); //10
58 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303

◆ reLinearizeTracks()

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 191 of file SequentialVertexFitter.cc.

192  {
193  GlobalPoint linP = state.position();
194  std::vector<RefCountedVertexTrack> finalTracks;
195  finalTracks.reserve(tracks.size());
196  for (typename std::vector<RefCountedVertexTrack>::const_iterator i = tracks.begin(); i != tracks.end(); i++) {
197  RefCountedLinearizedTrackState lTrData = (**i).linearizedTrack()->stateWithNewLinearizationPoint(linP);
198  // RefCountedLinearizedTrackState lTrData =
199  // theLTrackFactory->linearizedTrackState(linP,
200  // (**i).linearizedTrack()->track());
201  RefCountedVertexTrack vTrData = theVTrackFactory.vertexTrack(lTrData, state, (**i).weight());
202  finalTracks.push_back(vTrData);
203  }
204  return finalTracks;
205 }
VertexTrackFactory< N > theVTrackFactory
ReferenceCountingPointer< VertexTrack< N > > RefCountedVertexTrack
ReferenceCountingPointer< LinearizedTrackState< N > > RefCountedLinearizedTrackState

◆ setDefaultParameters()

template<unsigned int N>
void SequentialVertexFitter< N >::setDefaultParameters ( )
private

Definition at line 61 of file SequentialVertexFitter.cc.

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

61  {
62  thePSet.addParameter<double>("maxDistance", 0.01);
63  thePSet.addParameter<int>("maxNbrOfIterations", 10); //10
65 }
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:135

◆ setMaximumDistance()

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.

Referenced by GsfVertexFitter::GsfVertexFitter().

◆ setMaximumNumberOfIterations()

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.

Referenced by GsfVertexFitter::GsfVertexFitter().

◆ setTrackerBounds()

template<unsigned int N>
void SequentialVertexFitter< N >::setTrackerBounds ( float  radius,
float  halfLength 
)
inline

Definition at line 178 of file SequentialVertexFitter.h.

◆ vertex() [1/7]

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 68 of file SequentialVertexFitter.cc.

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

68  {
69  // Linearization Point
71  if (!insideTrackerBounds(linP))
72  linP = GlobalPoint(0, 0, 0);
73 
74  // Initial vertex state, with a very large error matrix
75  ROOT::Math::SMatrixIdentity id;
76  AlgebraicSymMatrix33 we(id);
77  GlobalError error(we * 10000);
78  VertexState state(linP, error);
79  std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
80  return fit(vtContainer, state, false);
81 }
const bool insideTrackerBounds(const GlobalPoint &point) const
LinearizationPointFinder * theLinP
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
virtual GlobalPoint getLinearizationPoint(const std::vector< reco::TransientTrack > &) const =0
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33

◆ vertex() [2/7]

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 91 of file SequentialVertexFitter.cc.

Referenced by Tau.Tau::dxy().

91  {
92  // Initial vertex state, with a very small weight matrix
93  GlobalPoint linP = tracks[0]->linearizedTrack()->linearizationPoint();
94  ROOT::Math::SMatrixIdentity id;
95  AlgebraicSymMatrix33 we(id);
96  GlobalError error(we * 10000);
97  VertexState state(linP, error);
98  return fit(tracks, state, false);
99 }
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33

◆ vertex() [3/7]

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 84 of file SequentialVertexFitter.cc.

Referenced by Tau.Tau::dxy().

85  {
86  VertexState state(spot);
87  return fit(tracks, state, true);
88 }
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const

◆ vertex() [4/7]

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 105 of file SequentialVertexFitter.cc.

Referenced by Tau.Tau::dxy().

106  {
107  // Initial vertex state, with a very large error matrix
108  ROOT::Math::SMatrixIdentity id;
109  AlgebraicSymMatrix33 we(id);
110  GlobalError error(we * 10000);
111  VertexState state(linPoint, error);
112  std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
113  return fit(vtContainer, state, false);
114 }
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33

◆ vertex() [5/7]

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 121 of file SequentialVertexFitter.cc.

Referenced by Tau.Tau::dxy().

122  {
123  VertexState beamSpotState(beamSpot);
124  std::vector<RefCountedVertexTrack> vtContainer;
125 
126  if (tracks.size() > 1) {
127  // Linearization Point search if there are more than 1 track
129  if (!insideTrackerBounds(linP))
130  linP = GlobalPoint(0, 0, 0);
131  ROOT::Math::SMatrixIdentity id;
132  AlgebraicSymMatrix33 we(id);
133  GlobalError error(we * 10000);
134  VertexState lpState(linP, error);
135  vtContainer = linearizeTracks(tracks, lpState);
136  } else {
137  // otherwise take the beamspot position.
138  vtContainer = linearizeTracks(tracks, beamSpotState);
139  }
140 
141  return fit(vtContainer, beamSpotState, true);
142 }
const bool insideTrackerBounds(const GlobalPoint &point) const
LinearizationPointFinder * theLinP
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
virtual GlobalPoint getLinearizationPoint(const std::vector< reco::TransientTrack > &) const =0
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33

◆ vertex() [6/7]

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 150 of file SequentialVertexFitter.cc.

Referenced by Tau.Tau::dxy().

152  {
153  VertexState state(priorPos, priorError);
154  std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
155  return fit(vtContainer, state, true);
156 }
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &tracks, const VertexState state) const

◆ vertex() [7/7]

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 163 of file SequentialVertexFitter.cc.

Referenced by Tau.Tau::dxy().

165  {
166  VertexState state(priorPos, priorError);
167  return fit(tracks, state, true);
168 }
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const

◆ vertexSmoother()

template<unsigned int N>
const VertexSmoother<N>* SequentialVertexFitter< N >::vertexSmoother ( ) const
inline

Definition at line 155 of file SequentialVertexFitter.h.

155 { return theSmoother; }
VertexSmoother< N > * theSmoother

◆ vertexUpdator()

template<unsigned int N>
const VertexUpdator<N>* SequentialVertexFitter< N >::vertexUpdator ( ) const
inline

Definition at line 153 of file SequentialVertexFitter.h.

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

Member Data Documentation

◆ theLinP

template<unsigned int N>
LinearizationPointFinder* SequentialVertexFitter< N >::theLinP
private

◆ theLTrackFactory

template<unsigned int N>
const AbstractLTSFactory<N>* SequentialVertexFitter< N >::theLTrackFactory
private

◆ theMaxShift

template<unsigned int N>
float SequentialVertexFitter< N >::theMaxShift
private

◆ theMaxStep

template<unsigned int N>
int SequentialVertexFitter< N >::theMaxStep
private

◆ thePSet

template<unsigned int N>
edm::ParameterSet SequentialVertexFitter< N >::thePSet
private

Definition at line 249 of file SequentialVertexFitter.h.

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

◆ theSmoother

template<unsigned int N>
VertexSmoother<N>* SequentialVertexFitter< N >::theSmoother
private

◆ theUpdator

template<unsigned int N>
VertexUpdator<N>* SequentialVertexFitter< N >::theUpdator
private

◆ theVTrackFactory

template<unsigned int N>
VertexTrackFactory<N> SequentialVertexFitter< N >::theVTrackFactory
private

Definition at line 255 of file SequentialVertexFitter.h.

◆ trackerBoundsHalfLength

template<unsigned int N>
float SequentialVertexFitter< N >::trackerBoundsHalfLength {273.5}
private

◆ trackerBoundsRadius

template<unsigned int N>
float SequentialVertexFitter< N >::trackerBoundsRadius {112.}
private