CMS 3D CMS Logo

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
 
 operator reco::VertexCompositePtrCandidate () 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
 
double priorTime () const
 
double priorTimeError () const
 
reco::TransientTrack refittedTrack (const reco::TransientTrack &track) const
 
std::vector< reco::TransientTrack > const & refittedTracks () const
 
void refittedTracks (const std::vector< reco::TransientTrack > &refittedTracks)
 
double time () const
 
double timeError () const
 
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 double time, 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 &pos, const double time, 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 double priorTime, const GlobalError &priorErr, const GlobalPoint &pos, const double time, 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 GlobalPoint &priorPos, const double priorTime, const GlobalError &priorErr, const GlobalPoint &pos, const double time, const GlobalError &posError, const std::vector< reco::TransientTrack > &tracks, float chi2, float ndf)
 
 TransientVertex (const GlobalPoint &priorPos, const GlobalError &priorErr, const double priorTime, const double priorTimeErr, const GlobalPoint &pos, const GlobalError &posError, const double time, const double timeErr, 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 18 of file TransientVertex.h.

Member Typedef Documentation

◆ TransientTrackToFloatMap

Definition at line 21 of file TransientVertex.h.

Constructor & Destructor Documentation

◆ TransientVertex() [1/14]

TransientVertex::TransientVertex ( )

Empty constructor, produces invalid vertex

Definition at line 11 of file TransientVertex.cc.

12  : theVertexState(),
14  theChi2(0),
15  theNDF(0),
16  vertexValid(false),
17  withPrior(false),
19  theCovMapAvailable(false),
20  withRefittedTracks(false) {}
std::vector< reco::TransientTrack > theOriginalTracks
VertexState theVertexState

◆ TransientVertex() [2/14]

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

References theNDF, and theOriginalTracks.

28  theChi2(chi2),
29  theNDF(0),
30  vertexValid(true),
31  withPrior(false),
33  theCovMapAvailable(false),
34  withRefittedTracks(false) {
35  theNDF = 2. * theOriginalTracks.size() - 3.;
36 }
std::vector< reco::TransientTrack > theOriginalTracks
auto const & tracks
cannot be loose
VertexState theVertexState

◆ TransientVertex() [3/14]

TransientVertex::TransientVertex ( const GlobalPoint pos,
const double  time,
const GlobalError posError,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2 
)

Definition at line 38 of file TransientVertex.cc.

References theNDF, and theOriginalTracks.

43  : theVertexState(pos, time, posTimeError),
45  theChi2(chi2),
46  theNDF(0),
47  vertexValid(true),
48  withPrior(false),
50  theCovMapAvailable(false),
51  withRefittedTracks(false) {
52  theNDF = 2. * theOriginalTracks.size() - 3.;
53 }
double time() const
std::vector< reco::TransientTrack > theOriginalTracks
auto const & tracks
cannot be loose
VertexState theVertexState

◆ TransientVertex() [4/14]

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 55 of file TransientVertex.cc.

62  theChi2(chi2),
63  theNDF(ndf),
64  vertexValid(true),
65  withPrior(false),
67  theCovMapAvailable(false),
68  withRefittedTracks(false) {}
std::vector< reco::TransientTrack > theOriginalTracks
auto const & tracks
cannot be loose
VertexState theVertexState

◆ TransientVertex() [5/14]

TransientVertex::TransientVertex ( const GlobalPoint pos,
const double  time,
const GlobalError posError,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2,
float  ndf 
)

Definition at line 70 of file TransientVertex.cc.

76  : theVertexState(pos, time, posTimeError),
78  theChi2(chi2),
79  theNDF(ndf),
80  vertexValid(true),
81  withPrior(false),
83  theCovMapAvailable(false),
84  withRefittedTracks(false) {}
double time() const
std::vector< reco::TransientTrack > theOriginalTracks
auto const & tracks
cannot be loose
VertexState theVertexState

◆ TransientVertex() [6/14]

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 86 of file TransientVertex.cc.

References theNDF, and theOriginalTracks.

92  : thePriorVertexState(priorPos, priorErr),
95  theChi2(chi2),
96  theNDF(0),
97  vertexValid(true),
98  withPrior(true),
100  theCovMapAvailable(false),
101  withRefittedTracks(false) {
102  theNDF = 2. * theOriginalTracks.size();
103 }
VertexState thePriorVertexState
std::vector< reco::TransientTrack > theOriginalTracks
auto const & tracks
cannot be loose
VertexState theVertexState

◆ TransientVertex() [7/14]

TransientVertex::TransientVertex ( const GlobalPoint priorPos,
const double  priorTime,
const GlobalError priorErr,
const GlobalPoint pos,
const double  time,
const GlobalError posError,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2 
)

Definition at line 105 of file TransientVertex.cc.

References theNDF, and theOriginalTracks.

113  : thePriorVertexState(priorPos, priorTime, priorErr),
116  theChi2(chi2),
117  theNDF(0),
118  vertexValid(true),
119  withPrior(true),
121  theCovMapAvailable(false),
122  withRefittedTracks(false) {
123  theNDF = 2. * theOriginalTracks.size();
124 }
double time() const
VertexState thePriorVertexState
std::vector< reco::TransientTrack > theOriginalTracks
double priorTime() const
auto const & tracks
cannot be loose
VertexState theVertexState

◆ TransientVertex() [8/14]

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 126 of file TransientVertex.cc.

133  : thePriorVertexState(priorPos, priorErr),
136  theChi2(chi2),
137  theNDF(ndf),
138  vertexValid(true),
139  withPrior(true),
141  theCovMapAvailable(false),
142  withRefittedTracks(false) {}
VertexState thePriorVertexState
std::vector< reco::TransientTrack > theOriginalTracks
auto const & tracks
cannot be loose
VertexState theVertexState

◆ TransientVertex() [9/14]

TransientVertex::TransientVertex ( const GlobalPoint priorPos,
const double  priorTime,
const GlobalError priorErr,
const GlobalPoint pos,
const double  time,
const GlobalError posError,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2,
float  ndf 
)

Definition at line 144 of file TransientVertex.cc.

153  : thePriorVertexState(priorPos, priorTime, priorErr),
156  theChi2(chi2),
157  theNDF(ndf),
158  vertexValid(true),
159  withPrior(true),
161  theCovMapAvailable(false),
162  withRefittedTracks(false) {}
double time() const
VertexState thePriorVertexState
std::vector< reco::TransientTrack > theOriginalTracks
double priorTime() const
auto const & tracks
cannot be loose
VertexState theVertexState

◆ TransientVertex() [10/14]

TransientVertex::TransientVertex ( const GlobalPoint priorPos,
const GlobalError priorErr,
const double  priorTime,
const double  priorTimeErr,
const GlobalPoint pos,
const GlobalError posError,
const double  time,
const double  timeErr,
const std::vector< reco::TransientTrack > &  tracks,
float  chi2,
float  ndf 
)

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

◆ TransientVertex() [11/14]

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 164 of file TransientVertex.cc.

References theNDF, and theOriginalTracks.

167  theChi2(chi2),
168  theNDF(0),
169  vertexValid(true),
170  withPrior(false),
172  theCovMapAvailable(false),
173  withRefittedTracks(false) {
174  theNDF = 2. * theOriginalTracks.size() - 3.;
175 }
std::vector< reco::TransientTrack > theOriginalTracks
auto const & tracks
cannot be loose
VertexState theVertexState

◆ TransientVertex() [12/14]

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 177 of file TransientVertex.cc.

183  theChi2(chi2),
184  theNDF(ndf),
185  vertexValid(true),
186  withPrior(false),
188  theCovMapAvailable(false),
189  withRefittedTracks(false) {}
std::vector< reco::TransientTrack > theOriginalTracks
auto const & tracks
cannot be loose
VertexState theVertexState

◆ TransientVertex() [13/14]

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

References theNDF, and theOriginalTracks.

198  theChi2(chi2),
199  theNDF(0),
200  vertexValid(true),
201  withPrior(true),
203  theCovMapAvailable(false),
204  withRefittedTracks(false) {
205  theNDF = 2. * theOriginalTracks.size();
206 }
VertexState thePriorVertexState
std::vector< reco::TransientTrack > theOriginalTracks
auto const & tracks
cannot be loose
VertexState theVertexState

◆ TransientVertex() [14/14]

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 208 of file TransientVertex.cc.

216  theChi2(chi2),
217  theNDF(ndf),
218  vertexValid(true),
219  withPrior(true),
221  theCovMapAvailable(false),
222  withRefittedTracks(false) {}
VertexState thePriorVertexState
std::vector< reco::TransientTrack > theOriginalTracks
auto const & tracks
cannot be loose
VertexState theVertexState

Member Function Documentation

◆ degreesOfFreedom()

float TransientVertex::degreesOfFreedom ( ) const
inline

◆ hasPrior()

bool TransientVertex::hasPrior ( ) const
inline

Definition at line 181 of file TransientVertex.h.

References withPrior.

Referenced by AdaptiveVertexReconstructor::cleanUp().

181 { return withPrior; }

◆ hasRefittedTracks()

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 206 of file TransientVertex.h.

References withRefittedTracks.

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

206 { return withRefittedTracks; }

◆ hasTrackWeight()

bool TransientVertex::hasTrackWeight ( ) const
inline

Returns true if the track-weights are available.

Definition at line 233 of file TransientVertex.h.

References theWeightMapIsAvailable.

233 { return theWeightMapIsAvailable; }

◆ isValid()

bool TransientVertex::isValid ( void  ) const
inline

◆ normalisedChiSquared()

float TransientVertex::normalisedChiSquared ( ) const
inline

Definition at line 189 of file TransientVertex.h.

References degreesOfFreedom(), and totalChiSquared().

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

189 { return totalChiSquared() / degreesOfFreedom(); }
float totalChiSquared() const
float degreesOfFreedom() const

◆ operator reco::Vertex()

TransientVertex::operator reco::Vertex ( ) const

Definition at line 304 of file TransientVertex.cc.

References mps_fire::i, sistrip::SpyUtilities::isValid(), trackWeight(), bphysicsOniaDQM_cfi::vertex, and HltBtagValidation_cff::Vertex.

304  {
305  //If the vertex is invalid, return an invalid TV !
306  if (!isValid())
307  return Vertex();
308 
310  // RecoVertex::convertError(theVertexState.error()),
313  totalChiSquared(),
315  theOriginalTracks.size());
316  for (std::vector<TransientTrack>::const_iterator i = theOriginalTracks.begin(); i != theOriginalTracks.end(); ++i) {
317  // const TrackTransientTrack* ttt = dynamic_cast<const TrackTransientTrack*>((*i).basicTransientTrack());
318  // if ((ttt!=0) && (ttt->persistentTrackRef().isNonnull()))
319  // {
320  // TrackRef tr = ttt->persistentTrackRef();
321  // TrackBaseRef tbr(tr);
322  if (withRefittedTracks) {
323  vertex.add((*i).trackBaseRef(), refittedTrack(*i).track(), trackWeight(*i));
324  } else {
325  vertex.add((*i).trackBaseRef(), trackWeight(*i));
326  }
327  //}
328  }
329  return vertex;
330 }
const Track & track() const
float totalChiSquared() const
reco::TransientTrack refittedTrack(const reco::TransientTrack &track) const
float degreesOfFreedom() const
const AlgebraicSymMatrix44 & matrix4D() const
bool isValid() const
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
std::vector< reco::TransientTrack > theOriginalTracks
double time() const
Definition: VertexState.h:73
VertexState theVertexState
float trackWeight(const reco::TransientTrack &track) const
GlobalError error4D() const
Definition: VertexState.h:67
GlobalPoint position() const
Definition: VertexState.h:62

◆ operator reco::VertexCompositePtrCandidate()

TransientVertex::operator reco::VertexCompositePtrCandidate ( ) const

Definition at line 332 of file TransientVertex.cc.

References reco::CompositePtrCandidate::addDaughter(), reco::CandidatePtrTransientTrack::candidate(), sistrip::SpyUtilities::isValid(), reco::Candidate::p4(), position, reco::VertexCompositePtrCandidate::setChi2AndNdof(), reco::VertexCompositePtrCandidate::setCovariance(), reco::LeafCandidate::setP4(), reco::VertexCompositePtrCandidate::setTime(), reco::LeafCandidate::setVertex(), protons_cff::time, trackWeight(), groupFilesInBlocks::tt, and x.

332  {
333  using namespace reco;
334  if (!isValid())
336 
337  VertexCompositePtrCandidate vtxCompPtrCand;
338 
339  vtxCompPtrCand.setTime(vertexState().time());
340  vtxCompPtrCand.setCovariance(vertexState().error4D().matrix4D());
341  vtxCompPtrCand.setChi2AndNdof(totalChiSquared(), degreesOfFreedom());
342  vtxCompPtrCand.setVertex(Candidate::Point(position().x(), position().y(), position().z()));
343 
345  for (std::vector<reco::TransientTrack>::const_iterator tt = theOriginalTracks.begin(); tt != theOriginalTracks.end();
346  ++tt) {
347  if (trackWeight(*tt) < 0.5)
348  continue;
349 
350  const CandidatePtrTransientTrack* cptt = dynamic_cast<const CandidatePtrTransientTrack*>(tt->basicTransientTrack());
351  if (cptt == nullptr)
352  edm::LogError("DynamicCastingFailed") << "Casting of TransientTrack to CandidatePtrTransientTrack failed!";
353  else {
354  p4 += cptt->candidate()->p4();
355  vtxCompPtrCand.addDaughter(cptt->candidate());
356  }
357  }
358 
359  //TODO: if has refitted tracks we should scale the candidate p4 to the refitted one
360  vtxCompPtrCand.setP4(p4);
361  return vtxCompPtrCand;
362 }
float totalChiSquared() const
double time() const
GlobalPoint position() const
CandidatePtr candidate() const override
Log< level::Error, false > LogError
float degreesOfFreedom() const
void setVertex(const Point &vertex) override
set vertex
VertexState const & vertexState() const
bool isValid() const
std::vector< reco::TransientTrack > theOriginalTracks
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
fixed size matrix
void setCovariance(const CovarianceMatrix &m)
set covariance matrix
void addDaughter(const CandidatePtr &)
add a daughter via a reference
float trackWeight(const reco::TransientTrack &track) const
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
void setChi2AndNdof(double chi2, double ndof)
set chi2 and ndof
void setP4(const LorentzVector &p4) final
set 4-momentum
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector

◆ originalTrack()

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 281 of file TransientVertex.cc.

References reco::TransientTrack::basicTransientTrack(), spr::find(), theOriginalTracks, and theRefittedTracks.

281  {
282  if (theRefittedTracks.empty())
283  throw VertexException("TransientVertex::requested No refitted tracks stored in vertex");
284  std::vector<TransientTrack>::const_iterator it = find(theRefittedTracks.begin(), theRefittedTracks.end(), refTrack);
285  if (it == theRefittedTracks.end())
286  throw VertexException(
287  "TransientVertex::requested Refitted track not found in list.\n address used for comparison: ")
288  << refTrack.basicTransientTrack();
289  size_t pos = it - theRefittedTracks.begin();
290  return theOriginalTracks[pos];
291 }
Common base class.
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::vector< reco::TransientTrack > theOriginalTracks
std::vector< reco::TransientTrack > theRefittedTracks
const BasicTransientTrack * basicTransientTrack() const

◆ originalTracks()

std::vector<reco::TransientTrack> const& TransientVertex::originalTracks ( ) const
inline

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

Definition at line 200 of file TransientVertex.h.

References theOriginalTracks.

Referenced by SplitVertexResolution::analyze(), AdaptiveVertexReconstructor::cleanUp(), AdaptiveVertexReconstructor::erase(), VertexFitterResult::fill(), VertexHigherPtSquared::operator()(), TrimmedVertexFitter::vertex(), and MultiVertexFitter::vertices().

200 { return theOriginalTracks; }
std::vector< reco::TransientTrack > theOriginalTracks

◆ position()

GlobalPoint TransientVertex::position ( ) const
inline

◆ positionError()

GlobalError TransientVertex::positionError ( ) const
inline

◆ priorError()

GlobalError TransientVertex::priorError ( ) const
inline

Definition at line 174 of file TransientVertex.h.

References VertexState::error(), VertexState::error4D(), VertexState::is4D(), and thePriorVertexState.

Referenced by AdaptiveVertexReconstructor::cleanUp().

174  {
176  }
VertexState thePriorVertexState
GlobalError error() const
Definition: VertexState.h:64
bool is4D() const
Definition: VertexState.h:95
GlobalError error4D() const
Definition: VertexState.h:67

◆ priorPosition()

GlobalPoint TransientVertex::priorPosition ( ) const
inline

Definition at line 173 of file TransientVertex.h.

References VertexState::position(), and thePriorVertexState.

Referenced by AdaptiveVertexReconstructor::cleanUp().

173 { return thePriorVertexState.position(); }
VertexState thePriorVertexState
GlobalPoint position() const
Definition: VertexState.h:62

◆ priorTime()

double TransientVertex::priorTime ( ) const
inline

Definition at line 179 of file TransientVertex.h.

References thePriorVertexState, and VertexState::time().

179 { return thePriorVertexState.time(); }
VertexState thePriorVertexState
double time() const
Definition: VertexState.h:73

◆ priorTimeError()

double TransientVertex::priorTimeError ( ) const
inline

Definition at line 180 of file TransientVertex.h.

References thePriorVertexState, and VertexState::timeError().

180 { return thePriorVertexState.timeError(); }
VertexState thePriorVertexState
double timeError() const
Definition: VertexState.h:75

◆ refittedTrack()

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 293 of file TransientVertex.cc.

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

Referenced by PFDisplacedVertexFinder::fitVertexFromSeed().

293  {
294  if (theRefittedTracks.empty())
295  throw VertexException("TransientVertex::requested No refitted tracks stored in vertex");
296  std::vector<TransientTrack>::const_iterator it = find(theOriginalTracks.begin(), theOriginalTracks.end(), track);
297  if (it == theOriginalTracks.end())
298  throw VertexException("transientVertex::requested Track not found in list.\n address used for comparison: ")
299  << track.basicTransientTrack();
300  size_t pos = it - theOriginalTracks.begin();
301  return theRefittedTracks[pos];
302 }
Common base class.
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::vector< reco::TransientTrack > theOriginalTracks
std::vector< reco::TransientTrack > theRefittedTracks

◆ refittedTracks() [1/2]

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 211 of file TransientVertex.h.

References theRefittedTracks.

Referenced by AdaptiveVertexReconstructor::cleanUp(), V0Fitter::fitAll(), CachingVertex< 5 >::operator TransientVertex(), PFTauSecondaryVertexProducer::produce(), and refittedTracks().

211 { return theRefittedTracks; }
std::vector< reco::TransientTrack > theRefittedTracks

◆ refittedTracks() [2/2]

void TransientVertex::refittedTracks ( const std::vector< reco::TransientTrack > &  refittedTracks)

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

Definition at line 229 of file TransientVertex.cc.

References refittedTracks(), theRefittedTracks, and withRefittedTracks.

229  {
230  if (refittedTracks.empty())
231  throw VertexException("TransientVertex::refittedTracks: No refitted tracks stored in input container");
233  withRefittedTracks = true;
234 }
Common base class.
std::vector< reco::TransientTrack > theRefittedTracks
std::vector< reco::TransientTrack > const & refittedTracks() const

◆ time()

double TransientVertex::time ( ) const
inline

Definition at line 177 of file TransientVertex.h.

References theVertexState, and VertexState::time().

177 { return theVertexState.time(); }
double time() const
Definition: VertexState.h:73
VertexState theVertexState

◆ timeError()

double TransientVertex::timeError ( ) const
inline

Definition at line 178 of file TransientVertex.h.

References theVertexState, and VertexState::timeError().

178 { return theVertexState.timeError(); }
double timeError() const
Definition: VertexState.h:75
VertexState theVertexState

◆ tkToTkCovariance() [1/2]

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 254 of file TransientVertex.cc.

References RandomServiceHelper::t1, RandomServiceHelper::t2, theCovMap, and theCovMapAvailable.

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

254  {
255  if (!theCovMapAvailable) {
256  throw VertexException("TransientVertex::Track-to-track covariance matrices not available");
257  }
258  const TransientTrack* tr1;
259  const TransientTrack* tr2;
260  if (t1 < t2) {
261  tr1 = &t1;
262  tr2 = &t2;
263  } else {
264  tr1 = &t2;
265  tr2 = &t1;
266  }
267  TTtoTTmap::const_iterator it = theCovMap.find(*tr1);
268  if (it != theCovMap.end()) {
269  const TTmap& tm = it->second;
270  TTmap::const_iterator nit = tm.find(*tr2);
271  if (nit != tm.end()) {
272  return (nit->second);
273  } else {
274  throw VertexException("TransientVertex::requested Track-to-track covariance matrix does not exist");
275  }
276  } else {
277  throw VertexException("TransientVertex::requested Track-to-track covariance matrix does not exist");
278  }
279 }
Common base class.
std::map< reco::TransientTrack, AlgebraicMatrix33 > TTmap
Definition: TTtoTTmap.h:11

◆ tkToTkCovariance() [2/2]

void TransientVertex::tkToTkCovariance ( const TTtoTTmap covMap)

Definition at line 236 of file TransientVertex.cc.

References theCovMap, and theCovMapAvailable.

236  {
237  theCovMap = covMap;
238  theCovMapAvailable = true;
239 }

◆ tkToTkCovarianceIsAvailable()

bool TransientVertex::tkToTkCovarianceIsAvailable ( ) const
inline

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

Definition at line 251 of file TransientVertex.h.

References theCovMapAvailable.

251 { return theCovMapAvailable; }

◆ totalChiSquared()

float TransientVertex::totalChiSquared ( ) const
inline

◆ trackWeight()

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 241 of file TransientVertex.cc.

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

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

241  {
243  std::vector<TransientTrack>::const_iterator foundTrack =
245  return ((foundTrack != theOriginalTracks.end()) ? 1. : 0.);
246  }
247  TransientTrackToFloatMap::const_iterator it = theWeightMap.find(track);
248  if (it != theWeightMap.end()) {
249  return (it->second);
250  }
251  return 0.;
252 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
TransientTrackToFloatMap theWeightMap
std::vector< reco::TransientTrack > theOriginalTracks

◆ vertexState()

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 168 of file TransientVertex.h.

References theVertexState.

Referenced by AdaptiveVertexReconstructor::cleanUp().

168 { return theVertexState; }
VertexState theVertexState

◆ weightMap() [1/2]

TransientTrackToFloatMap TransientVertex::weightMap ( ) const
inline

Definition at line 244 of file TransientVertex.h.

References theWeightMap.

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

244 { return theWeightMap; }
TransientTrackToFloatMap theWeightMap

◆ weightMap() [2/2]

void TransientVertex::weightMap ( const TransientTrackToFloatMap theMap)

Definition at line 224 of file TransientVertex.cc.

References theWeightMap, and theWeightMapIsAvailable.

224  {
225  theWeightMap = theMap;
227 }
TransientTrackToFloatMap theWeightMap

Member Data Documentation

◆ theChi2

float TransientVertex::theChi2
private

Definition at line 273 of file TransientVertex.h.

Referenced by totalChiSquared().

◆ theCovMap

TTtoTTmap TransientVertex::theCovMap
private

Definition at line 278 of file TransientVertex.h.

Referenced by tkToTkCovariance().

◆ theCovMapAvailable

bool TransientVertex::theCovMapAvailable
private

Definition at line 276 of file TransientVertex.h.

Referenced by tkToTkCovariance(), and tkToTkCovarianceIsAvailable().

◆ theNDF

float TransientVertex::theNDF
private

Definition at line 274 of file TransientVertex.h.

Referenced by degreesOfFreedom(), and TransientVertex().

◆ theOriginalTracks

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

◆ thePriorVertexState

VertexState TransientVertex::thePriorVertexState
mutableprivate

Definition at line 265 of file TransientVertex.h.

Referenced by priorError(), priorPosition(), priorTime(), and priorTimeError().

◆ theRefittedTracks

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

Definition at line 271 of file TransientVertex.h.

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

◆ theVertexState

VertexState TransientVertex::theVertexState
mutableprivate

Definition at line 266 of file TransientVertex.h.

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

◆ theWeightMap

TransientTrackToFloatMap TransientVertex::theWeightMap
private

Definition at line 279 of file TransientVertex.h.

Referenced by trackWeight(), and weightMap().

◆ theWeightMapIsAvailable

bool TransientVertex::theWeightMapIsAvailable
private

Definition at line 276 of file TransientVertex.h.

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

◆ vertexValid

bool TransientVertex::vertexValid
private

Definition at line 275 of file TransientVertex.h.

Referenced by isValid().

◆ withPrior

bool TransientVertex::withPrior
private

Definition at line 276 of file TransientVertex.h.

Referenced by hasPrior().

◆ withRefittedTracks

bool TransientVertex::withRefittedTracks
private

Definition at line 277 of file TransientVertex.h.

Referenced by hasRefittedTracks(), and refittedTracks().