CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Private Attributes
TransientVertex Class Reference

#include <TransientVertex.h>

Public Types

typedef std::map
< reco::TransientTrack, float > 
TransientTrackToFloatMap
 

Public Member Functions

float degreesOfFreedom () const
 
bool hasPrior () const
 
bool hasRefittedTracks () const
 
bool hasTrackWeight () const
 
bool isValid () const
 
float normalisedChiSquared () const
 
 operator reco::Vertex () const
 
reco::TransientTrack originalTrack (const reco::TransientTrack &refTrack) const
 
std::vector
< reco::TransientTrack > const & 
originalTracks () const
 
GlobalPoint position () const
 
GlobalError positionError () const
 
GlobalError priorError () const
 
GlobalPoint priorPosition () const
 
reco::TransientTrack refittedTrack (const reco::TransientTrack &track) const
 
std::vector
< reco::TransientTrack > const & 
refittedTracks () const
 
void refittedTracks (const std::vector< reco::TransientTrack > &refittedTracks)
 
AlgebraicMatrix33 tkToTkCovariance (const reco::TransientTrack &t1, const reco::TransientTrack &t2) const
 
void tkToTkCovariance (const TTtoTTmap &covMap)
 
bool tkToTkCovarianceIsAvailable () const
 
float totalChiSquared () const
 
float trackWeight (const reco::TransientTrack &track) const
 
 TransientVertex ()
 
 TransientVertex (const GlobalPoint &pos, const GlobalError &posError, const std::vector< reco::TransientTrack > &tracks, float chi2)
 
 TransientVertex (const GlobalPoint &pos, const GlobalError &posError, const std::vector< reco::TransientTrack > &tracks, float chi2, float ndf)
 
 TransientVertex (const GlobalPoint &priorPos, const GlobalError &priorErr, const GlobalPoint &pos, const GlobalError &posError, const std::vector< reco::TransientTrack > &tracks, float chi2)
 
 TransientVertex (const GlobalPoint &priorPos, const GlobalError &priorErr, const GlobalPoint &pos, const GlobalError &posError, const std::vector< reco::TransientTrack > &tracks, float chi2, float ndf)
 
 TransientVertex (const VertexState &state, const std::vector< reco::TransientTrack > &tracks, float chi2)
 
 TransientVertex (const VertexState &state, const std::vector< reco::TransientTrack > &tracks, float chi2, float ndf)
 
 TransientVertex (const VertexState &prior, const VertexState &state, const std::vector< reco::TransientTrack > &tracks, float chi2)
 
 TransientVertex (const VertexState &prior, const VertexState &state, const std::vector< reco::TransientTrack > &tracks, float chi2, float ndf)
 
VertexState const & vertexState () const
 
TransientTrackToFloatMap weightMap () const
 
void weightMap (const TransientTrackToFloatMap &theMap)
 

Private Attributes

float theChi2
 
TTtoTTmap theCovMap
 
bool theCovMapAvailable
 
float theNDF
 
std::vector< reco::TransientTracktheOriginalTracks
 
VertexState thePriorVertexState
 
std::vector< reco::TransientTracktheRefittedTracks
 
VertexState theVertexState
 
TransientTrackToFloatMap theWeightMap
 
bool theWeightMapIsAvailable
 
bool vertexValid
 
bool withPrior
 
bool withRefittedTracks
 

Detailed Description

Definition at line 17 of file TransientVertex.h.

Member Typedef Documentation

Definition at line 21 of file TransientVertex.h.

Constructor & Destructor Documentation

TransientVertex::TransientVertex ( )

Empty constructor, produces invalid vertex

Definition at line 9 of file TransientVertex.cc.

10  theChi2(0), theNDF(0), vertexValid(false), withPrior(false),
12  withRefittedTracks(false)
13 {}
std::vector< reco::TransientTrack > theOriginalTracks
VertexState theVertexState
TransientVertex::TransientVertex ( const GlobalPoint pos,
const GlobalError posError,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2 
)

Constructor defining the RecVertex by its 3D position and position uncertainty, its associated tracks and its chi-squared. The number of degrees of freedom is equal to 2*nb of tracks - 3.

Definition at line 16 of file TransientVertex.cc.

References theNDF, and theOriginalTracks.

17  :
18  theVertexState(pos, posError), theOriginalTracks(tracks),
19  theChi2(chi2), theNDF(0), vertexValid(true), withPrior(false),
21  withRefittedTracks(false)
22 {
23  theNDF = 2.*theOriginalTracks.size() - 3.;
24  // addTracks(tracks);
25 }
std::vector< reco::TransientTrack > theOriginalTracks
tuple tracks
Definition: testEve_cfg.py:39
VertexState theVertexState
TransientVertex::TransientVertex ( const GlobalPoint pos,
const GlobalError posError,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2,
float  ndf 
)

Constructor defining the RecVertex by its 3D position and position uncertainty, its associated tracks, its chi-squared and its number of degrees of freedom. The ndf can be a float.

Definition at line 28 of file TransientVertex.cc.

29  :
30  theVertexState(pos, posError), theOriginalTracks(tracks),
31  theChi2(chi2), theNDF(ndf), vertexValid(true), withPrior(false),
33  withRefittedTracks(false)
34 {
35  // addTracks(tracks);
36 }
std::vector< reco::TransientTrack > theOriginalTracks
tuple tracks
Definition: testEve_cfg.py:39
VertexState theVertexState
TransientVertex::TransientVertex ( const GlobalPoint priorPos,
const GlobalError priorErr,
const GlobalPoint pos,
const GlobalError posError,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2 
)

Constructor defining the RecVertex by the prior, the vertex 3D position and uncertainty, the associated tracks and the chi-squared. Since the prior brings information on 3 coordinates, the number of degrees of freedom is equal to 2*nb of tracks.

Definition at line 39 of file TransientVertex.cc.

References theNDF, and theOriginalTracks.

41  :
42  thePriorVertexState(priorPos, priorErr), theVertexState(pos, posError),
45  withRefittedTracks(false)
46 {
47  theNDF = 2.*theOriginalTracks.size();
48  // addTracks(tracks);
49 }
VertexState thePriorVertexState
std::vector< reco::TransientTrack > theOriginalTracks
tuple tracks
Definition: testEve_cfg.py:39
VertexState theVertexState
TransientVertex::TransientVertex ( const GlobalPoint priorPos,
const GlobalError priorErr,
const GlobalPoint pos,
const GlobalError posError,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2,
float  ndf 
)

Constructor defining the RecVertex by the prior, the vertex 3D position and uncertainty, the associated tracks, the chi-squared and the number of degrees of freedom. The ndf can be a float.

Definition at line 52 of file TransientVertex.cc.

54  :
55  thePriorVertexState(priorPos, priorErr), theVertexState(pos, posError),
56  theOriginalTracks(tracks), theChi2(chi2), theNDF(ndf), vertexValid(true),
58  withRefittedTracks(false)
59 {
60  // addTracks(tracks);
61 }
VertexState thePriorVertexState
std::vector< reco::TransientTrack > theOriginalTracks
tuple tracks
Definition: testEve_cfg.py:39
VertexState theVertexState
TransientVertex::TransientVertex ( const VertexState &  state,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2 
)

Constructor defining the RecVertex by its 3D position and position uncertainty, its associated tracks and its chi-squared. The number of degrees of freedom is equal to 2*nb of tracks - 3.

Definition at line 64 of file TransientVertex.cc.

References theNDF, and theOriginalTracks.

65  :
67  theChi2(chi2), theNDF(0), vertexValid(true), withPrior(false),
69  withRefittedTracks(false)
70 {
71  theNDF = 2.*theOriginalTracks.size() - 3.;
72 }
std::vector< reco::TransientTrack > theOriginalTracks
tuple tracks
Definition: testEve_cfg.py:39
VertexState theVertexState
TransientVertex::TransientVertex ( const VertexState &  state,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2,
float  ndf 
)

Constructor defining the RecVertex by its 3D position and position uncertainty, its associated tracks, its chi-squared and its number of degrees of freedom. The ndf can be a float.

Definition at line 75 of file TransientVertex.cc.

76  :
78  theChi2(chi2), theNDF(ndf), vertexValid(true), withPrior(false),
80  withRefittedTracks(false)
81 {
82  // addTracks(tracks);
83 }
std::vector< reco::TransientTrack > theOriginalTracks
tuple tracks
Definition: testEve_cfg.py:39
VertexState theVertexState
TransientVertex::TransientVertex ( const VertexState &  prior,
const VertexState &  state,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2 
)

Constructor defining the RecVertex by the prior, the vertex 3D position and uncertainty, the associated tracks and the chi-squared. Since the prior brings information on 3 coordinates, the number of degrees of freedom is equal to 2*nb of tracks.

Definition at line 86 of file TransientVertex.cc.

References theNDF, and theOriginalTracks.

89  :
93  withRefittedTracks(false)
94 {
95  theNDF = 2.*theOriginalTracks.size();
96  // addTracks(tracks);
97 }
VertexState thePriorVertexState
std::vector< reco::TransientTrack > theOriginalTracks
tuple tracks
Definition: testEve_cfg.py:39
dictionary prior
VertexState theVertexState
TransientVertex::TransientVertex ( const VertexState &  prior,
const VertexState &  state,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2,
float  ndf 
)

Constructor defining the RecVertex by the prior, the vertex 3D position and uncertainty, the associated tracks, the chi-squared and the number of degrees of freedom. The ndf can be a float.

Definition at line 100 of file TransientVertex.cc.

103  :
105  theOriginalTracks(tracks), theChi2(chi2), theNDF(ndf), vertexValid(true),
106  withPrior(true), theWeightMapIsAvailable(false),
108 {
109  // addTracks(tracks);
110 }
VertexState thePriorVertexState
std::vector< reco::TransientTrack > theOriginalTracks
tuple tracks
Definition: testEve_cfg.py:39
dictionary prior
VertexState theVertexState

Member Function Documentation

float TransientVertex::degreesOfFreedom ( ) const
inline
bool TransientVertex::hasPrior ( ) const
inline

Definition at line 117 of file TransientVertex.h.

References withPrior.

Referenced by AdaptiveVertexReconstructor::cleanUp().

117 { return withPrior; }
bool TransientVertex::hasRefittedTracks ( ) const
inline

Returns true if at for at least one of the original tracks, the refitted track is available

Definition at line 149 of file TransientVertex.h.

References withRefittedTracks.

Referenced by AdaptiveVertexReconstructor::cleanUp(), V0Fitter::fitAll(), PFTauSecondaryVertexProducer::produce(), and HPSPFRecoTauAlgorithm::refitThreeProng().

149 { return withRefittedTracks; }
bool TransientVertex::hasTrackWeight ( ) const
inline

Returns true if the track-weights are available.

Definition at line 180 of file TransientVertex.h.

References theWeightMapIsAvailable.

Referenced by PrimaryVertexValidation::analyze().

180 { return theWeightMapIsAvailable; }
bool TransientVertex::isValid ( void  ) const
inline
float TransientVertex::normalisedChiSquared ( ) const
inline

Definition at line 125 of file TransientVertex.h.

References degreesOfFreedom(), and totalChiSquared().

Referenced by PrimaryVertexValidation::analyze(), HLTmmkkFilter::hltFilter(), and HLTmmkFilter::hltFilter().

125  {
126  return totalChiSquared() / degreesOfFreedom();
127  }
float totalChiSquared() const
float degreesOfFreedom() const
TransientVertex::operator reco::Vertex ( ) const

Definition at line 206 of file TransientVertex.cc.

References reco::Vertex::add(), and i.

207 {
208  //If the vertex is invalid, return an invalid TV !
209  if (!isValid()) return Vertex();
210 
211  Vertex vertex(Vertex::Point(theVertexState.position()),
212 // RecoVertex::convertError(theVertexState.error()),
213  theVertexState.error().matrix_new(),
215  for (std::vector<TransientTrack>::const_iterator i = theOriginalTracks.begin();
216  i != theOriginalTracks.end(); ++i) {
217 // const TrackTransientTrack* ttt = dynamic_cast<const TrackTransientTrack*>((*i).basicTransientTrack());
218 // if ((ttt!=0) && (ttt->persistentTrackRef().isNonnull()))
219 // {
220 // TrackRef tr = ttt->persistentTrackRef();
221 // TrackBaseRef tbr(tr);
222  if (withRefittedTracks) {
223 
224  vertex.add((*i).trackBaseRef(), refittedTrack(*i).track(), trackWeight ( *i ) );
225  } else {
226  vertex.add((*i).trackBaseRef(), trackWeight ( *i ) );
227  }
228  //}
229  }
230  return vertex;
231 }
int i
Definition: DBlmapReader.cc:9
float totalChiSquared() const
reco::TransientTrack refittedTrack(const reco::TransientTrack &track) const
float degreesOfFreedom() const
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
std::vector< reco::TransientTrack > theOriginalTracks
float trackWeight(const reco::TransientTrack &track) const
void add(const TrackBaseRef &r, float w=1.0)
add a reference to a Track
const Track & track() const
VertexState theVertexState
bool isValid() const
TransientTrack TransientVertex::originalTrack ( const reco::TransientTrack refTrack) const

Returns the original track which corresponds to a particular refitted Track Throws an exception if now refitted tracks are stored ot the track is not found in the list

Definition at line 182 of file TransientVertex.cc.

References spr::find(), theOriginalTracks, and theRefittedTracks.

183 {
184  if (theRefittedTracks.empty())
185  throw VertexException("No refitted tracks stored in vertex");
186  std::vector<TransientTrack>::const_iterator it =
187  find(theRefittedTracks.begin(), theRefittedTracks.end(), refTrack);
188  if (it==theRefittedTracks.end())
189  throw VertexException("Refitted track not found in list");
190  size_t pos = it - theRefittedTracks.begin();
191  return theOriginalTracks[pos];
192 }
Common base class.
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
std::vector< reco::TransientTrack > theOriginalTracks
std::vector< reco::TransientTrack > theRefittedTracks
std::vector<reco::TransientTrack> const& TransientVertex::originalTracks ( ) const
inline
GlobalPoint TransientVertex::position ( ) const
inline
GlobalError TransientVertex::positionError ( ) const
inline
GlobalError TransientVertex::priorError ( ) const
inline

Definition at line 116 of file TransientVertex.h.

References thePriorVertexState.

Referenced by AdaptiveVertexReconstructor::cleanUp().

116 { return thePriorVertexState.error(); }
VertexState thePriorVertexState
GlobalPoint TransientVertex::priorPosition ( ) const
inline

Definition at line 115 of file TransientVertex.h.

References thePriorVertexState.

Referenced by AdaptiveVertexReconstructor::cleanUp().

115 { return thePriorVertexState.position(); }
VertexState thePriorVertexState
TransientTrack TransientVertex::refittedTrack ( const reco::TransientTrack track) const

Returns the refitted track which corresponds to a particular original Track Throws an exception if now refitted tracks are stored ot the track is not found in the list

Definition at line 194 of file TransientVertex.cc.

References spr::find(), theOriginalTracks, and theRefittedTracks.

Referenced by PFDisplacedVertexFinder::fitVertexFromSeed().

195 {
196  if (theRefittedTracks.empty())
197  throw VertexException("No refitted tracks stored in vertex");
198  std::vector<TransientTrack>::const_iterator it =
199  find(theOriginalTracks.begin(), theOriginalTracks.end(), track);
200  if (it==theOriginalTracks.end())
201  throw VertexException("Track not found in list");
202  size_t pos = it - theOriginalTracks.begin();
203  return theRefittedTracks[pos];
204 }
Common base class.
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
std::vector< reco::TransientTrack > theOriginalTracks
std::vector< reco::TransientTrack > theRefittedTracks
std::vector<reco::TransientTrack> const& TransientVertex::refittedTracks ( ) const
inline

Access to the refitted tracks used to make the vertex. Returns track container by value.

Definition at line 155 of file TransientVertex.h.

References theRefittedTracks.

Referenced by AdaptiveVertexReconstructor::cleanUp(), V0Fitter::fitAll(), CachingVertex< N >::operator TransientVertex(), PFTauSecondaryVertexProducer::produce(), refittedTracks(), HPSPFRecoTauAlgorithm::refitThreeProng(), CandCommonVertexFitterBase::set(), and PFCandCommonVertexFitterBase::set().

155  {
156  return theRefittedTracks;
157  }
std::vector< reco::TransientTrack > theRefittedTracks
void TransientVertex::refittedTracks ( const std::vector< reco::TransientTrack > &  refittedTracks)

Method to set the refitted tracks used to make the vertex.

Definition at line 120 of file TransientVertex.cc.

References refittedTracks(), theRefittedTracks, and withRefittedTracks.

122 {
123  if (refittedTracks.empty())
124  throw VertexException("TransientVertex::refittedTracks: No refitted tracks stored in input container");
126  withRefittedTracks = true;
127 }
Common base class.
std::vector< reco::TransientTrack > const & refittedTracks() const
std::vector< reco::TransientTrack > theRefittedTracks
AlgebraicMatrix33 TransientVertex::tkToTkCovariance ( const reco::TransientTrack t1,
const reco::TransientTrack t2 
) const

Returns the Track-to-track covariance matrix for two specified tracks. In case these do not exist, or one of the tracks does not belong to the vertex, an exception is thrown.

Definition at line 151 of file TransientVertex.cc.

References reco::t2, theCovMap, and theCovMapAvailable.

Referenced by CachingVertex< N >::operator TransientVertex().

152 {
153  if (!theCovMapAvailable) {
154  throw VertexException("TransientVertex::Track-to-track covariance matrices not available");
155  }
156  const TransientTrack* tr1;
157  const TransientTrack* tr2;
158  if (t1<t2) {
159  tr1 = &t1;
160  tr2 = &t2;
161  } else {
162  tr1 = &t2;
163  tr2 = &t1;
164  }
165  TTtoTTmap::const_iterator it = theCovMap.find(*tr1);
166  if (it != theCovMap.end()) {
167  const TTmap & tm = it->second;
168  TTmap::const_iterator nit = tm.find(*tr2);
169  if (nit != tm.end()) {
170  return( nit->second);
171  }
172  else {
173  throw VertexException("TransientVertex::requested Track-to-track covariance matrix does not exist");
174  }
175  }
176  else {
177  throw VertexException("TransientVertex::requested Track-to-track covariance matrix does not exist");
178  }
179 }
Common base class.
std::map< reco::TransientTrack, AlgebraicMatrix33 > TTmap
Definition: TTtoTTmap.h:11
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:18
void TransientVertex::tkToTkCovariance ( const TTtoTTmap covMap)

Definition at line 130 of file TransientVertex.cc.

References theCovMap, and theCovMapAvailable.

131 {
132  theCovMap = covMap;
133  theCovMapAvailable = true;
134 }
bool TransientVertex::tkToTkCovarianceIsAvailable ( ) const
inline

Returns true if the Track-to-track covariance matrices have been calculated.

Definition at line 199 of file TransientVertex.h.

References theCovMapAvailable.

199 { return theCovMapAvailable; }
float TransientVertex::totalChiSquared ( ) const
inline

Implements method of abstract Vertex. Returns track pointer container by value

Definition at line 124 of file TransientVertex.h.

References theChi2.

Referenced by PrimaryVertexValidation::analyze(), AdaptiveVertexReconstructor::cleanUp(), VertexFitterResult::fill(), V0Fitter::fitAll(), PFDisplacedVertexFinder::fitVertexFromSeed(), normalisedChiSquared(), and TrimmedVertexFitter::vertex().

124 { return theChi2; }
float TransientVertex::trackWeight ( const reco::TransientTrack track) const

Returns the weight with which a track has been used in the fit. If the track is not present in the list, it is obviously not used, and this method returns a weight of 0. If this information has not been provided at construction, a weight of 1.0 is assumed for all tracks used and present in the originalTracks() std::vector.

Definition at line 136 of file TransientVertex.cc.

References spr::find(), theOriginalTracks, theWeightMap, and theWeightMapIsAvailable.

Referenced by PrimaryVertexValidation::analyze(), AdaptiveVertexReconstructor::cleanUp(), AdaptiveVertexReconstructor::erase(), and PFDisplacedVertexFinder::fitVertexFromSeed().

136  {
138  std::vector<TransientTrack>::const_iterator foundTrack = find(theOriginalTracks.begin(),
139  theOriginalTracks.end(), track);
140  return ((foundTrack != theOriginalTracks.end()) ? 1. : 0.);
141  }
142  TransientTrackToFloatMap::const_iterator it = theWeightMap.find(track);
143  if (it != theWeightMap.end()) {
144  return(it->second);
145  }
146  return 0.;
147 
148 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
TransientTrackToFloatMap theWeightMap
std::vector< reco::TransientTrack > theOriginalTracks
VertexState const& TransientVertex::vertexState ( ) const
inline

Constructor defining the RecVertex by its 3D position and position uncertainty, its associated tracks, its chi-squared and its number of degrees of freedom, and the track weights. The ndf can be a float.Access methods

Definition at line 112 of file TransientVertex.h.

References theVertexState.

Referenced by AdaptiveVertexReconstructor::cleanUp().

112 { return theVertexState; }
VertexState theVertexState
TransientTrackToFloatMap TransientVertex::weightMap ( ) const
inline
void TransientVertex::weightMap ( const TransientTrackToFloatMap theMap)

Definition at line 112 of file TransientVertex.cc.

References theWeightMap, and theWeightMapIsAvailable.

113 {
114  theWeightMap = theMap;
116 // removeTracks(); // remove trackrefs from reco::Vertex
117 // addTracks( theOriginalTracks );
118 }
TransientTrackToFloatMap theWeightMap

Member Data Documentation

float TransientVertex::theChi2
private

Definition at line 224 of file TransientVertex.h.

Referenced by totalChiSquared().

TTtoTTmap TransientVertex::theCovMap
private

Definition at line 229 of file TransientVertex.h.

Referenced by tkToTkCovariance().

bool TransientVertex::theCovMapAvailable
private

Definition at line 227 of file TransientVertex.h.

Referenced by tkToTkCovariance(), and tkToTkCovarianceIsAvailable().

float TransientVertex::theNDF
private

Definition at line 225 of file TransientVertex.h.

Referenced by degreesOfFreedom(), and TransientVertex().

std::vector<reco::TransientTrack> TransientVertex::theOriginalTracks
private
VertexState TransientVertex::thePriorVertexState
mutableprivate

Definition at line 215 of file TransientVertex.h.

Referenced by priorError(), and priorPosition().

std::vector<reco::TransientTrack> TransientVertex::theRefittedTracks
private

Definition at line 221 of file TransientVertex.h.

Referenced by originalTrack(), refittedTrack(), and refittedTracks().

VertexState TransientVertex::theVertexState
mutableprivate

Definition at line 216 of file TransientVertex.h.

Referenced by position(), positionError(), and vertexState().

TransientTrackToFloatMap TransientVertex::theWeightMap
private

Definition at line 230 of file TransientVertex.h.

Referenced by trackWeight(), and weightMap().

bool TransientVertex::theWeightMapIsAvailable
private

Definition at line 227 of file TransientVertex.h.

Referenced by hasTrackWeight(), trackWeight(), and weightMap().

bool TransientVertex::vertexValid
private

Definition at line 226 of file TransientVertex.h.

Referenced by isValid().

bool TransientVertex::withPrior
private

Definition at line 227 of file TransientVertex.h.

Referenced by hasPrior().

bool TransientVertex::withRefittedTracks
private

Definition at line 228 of file TransientVertex.h.

Referenced by hasRefittedTracks(), and refittedTracks().