CMS 3D CMS Logo

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

#include <Vertex.h>

Inheritance diagram for reco::Vertex:
reco::PFDisplacedVertex reco::SecondaryVertex

Classes

class  TrackEqual
 

Public Types

enum  { dimension = 3 }
 error matrix dimension More...
 
enum  { size = dimension * ( dimension + 1 ) / 2 }
 matix size More...
 
typedef math::Error< dimension >
::type 
CovarianceMatrix
 covariance error matrix (3x3) More...
 
typedef math::Error< dimension >
::type 
Error
 covariance error matrix (3x3) More...
 
typedef unsigned int index
 index type More...
 
typedef math::XYZPoint Point
 point in the space More...
 
typedef std::vector
< TrackBaseRef >
::const_iterator 
trackRef_iterator
 The iteratator for the vector<TrackRef> More...
 

Public Member Functions

void add (const TrackBaseRef &r, float w=1.0)
 add a reference to a Track More...
 
void add (const TrackBaseRef &r, const Track &refTrack, float w=1.0)
 add the original a Track(reference) and the smoothed Track More...
 
double chi2 () const
 chi-squares More...
 
double covariance (int i, int j) const
 (i, j)-th element of error matrix, i, j = 0, ... 2 More...
 
CovarianceMatrix covariance () const
 return SMatrix More...
 
Error error () const
 return SMatrix More...
 
void fill (CovarianceMatrix &v) const
 fill SMatrix More...
 
bool hasRefittedTracks () const
 Checks whether refitted tracks are stored. More...
 
bool isFake () const
 
bool isValid () const
 Tells whether the vertex is valid. More...
 
double ndof () const
 
double normalizedChi2 () const
 chi-squared divided by n.d.o.f. More...
 
unsigned int nTracks (float minWeight=0.5) const
 Returns the number of tracks in the vertex with weight above minWeight. More...
 
TrackBaseRef originalTrack (const Track &refTrack) const
 
math::XYZTLorentzVectorD p4 (float mass=0.13957018, float minWeight=0.5) const
 Returns the four momentum of the sum of the tracks, assuming the given mass for the decay products. More...
 
const Pointposition () const
 position More...
 
Track refittedTrack (const TrackBaseRef &track) const
 
Track refittedTrack (const TrackRef &track) const
 
const std::vector< Track > & refittedTracks () const
 Returns the container of refitted tracks. More...
 
void removeTracks ()
 
trackRef_iterator tracks_begin () const
 first iterator over tracks More...
 
trackRef_iterator tracks_end () const
 last iterator over tracks More...
 
size_t tracksSize () const
 number of tracks More...
 
float trackWeight (const TrackBaseRef &r) const
 returns the weight with which a Track has contributed to the vertex-fit. More...
 
float trackWeight (const TrackRef &r) const
 returns the weight with which a Track has contributed to the vertex-fit. More...
 
 Vertex ()
 
 Vertex (const Point &, const Error &)
 Constructor for a fake vertex. More...
 
 Vertex (const Point &, const Error &, double chi2, double ndof, size_t size)
 constructor for a valid vertex, with all data More...
 
double x () const
 x coordinate More...
 
double xError () const
 error on x More...
 
double y () const
 y coordinate More...
 
double yError () const
 error on y More...
 
double z () const
 y coordinate More...
 
double zError () const
 error on z More...
 

Private Member Functions

index idx (index i, index j) const
 position index More...
 

Private Attributes

float chi2_
 chi-sqared More...
 
float covariance_ [size]
 covariance matrix (3x3) as vector More...
 
float ndof_
 number of degrees of freedom More...
 
Point position_
 position More...
 
std::vector< TrackrefittedTracks_
 The vector of refitted tracks. More...
 
std::vector< TrackBaseReftracks_
 reference to tracks More...
 
bool validity_
 tells wether the vertex is really valid. More...
 
std::vector< uint8_t > weights_
 

Detailed Description

A reconstructed Vertex providing position, error, chi2, ndof and reconstrudted tracks. The vertex can be valid, fake, or invalid. A valid vertex is one which has been obtained from a vertex fit of tracks, and all data is meaningful A fake vertex is a vertex which was not made out of a proper fit with tracks, but still has a position and error (chi2 and ndof are null). For a primary vertex, it could simply be the beam line. A fake vertex is considered valid. An invalid vertex has no meaningful data.

Author
Luca Lista, INFN
Version
Id:
Vertex.h,v 1.39 2011/02/25 22:05:02 dlange Exp

Definition at line 35 of file Vertex.h.

Member Typedef Documentation

covariance error matrix (3x3)

Definition at line 46 of file Vertex.h.

covariance error matrix (3x3)

Definition at line 44 of file Vertex.h.

typedef unsigned int reco::Vertex::index

index type

Definition at line 50 of file Vertex.h.

point in the space

Definition at line 40 of file Vertex.h.

typedef std::vector<TrackBaseRef >::const_iterator reco::Vertex::trackRef_iterator

The iteratator for the vector<TrackRef>

Definition at line 38 of file Vertex.h.

Member Enumeration Documentation

anonymous enum

error matrix dimension

Enumerator
dimension 

Definition at line 42 of file Vertex.h.

anonymous enum

matix size

Enumerator
size 

Definition at line 48 of file Vertex.h.

Constructor & Destructor Documentation

reco::Vertex::Vertex ( )
inline

default constructor - The vertex will not be valid. Position, error, chi2, ndof will have random entries, and the vectors of tracks will be empty Use the isValid method to check that your vertex is valid.

Definition at line 54 of file Vertex.h.

References covariance_, i, size, and validity_.

54  : chi2_( 0.0 ), ndof_( 0 ), position_(0.,0.,0. ) { validity_ = false; for(int i = 0; i < size; ++ i ) covariance_[i]=0.;
55 }
int i
Definition: DBlmapReader.cc:9
float chi2_
chi-sqared
Definition: Vertex.h:153
Point position_
position
Definition: Vertex.h:157
float covariance_[size]
covariance matrix (3x3) as vector
Definition: Vertex.h:159
float ndof_
number of degrees of freedom
Definition: Vertex.h:155
bool validity_
tells wether the vertex is really valid.
Definition: Vertex.h:166
Vertex::Vertex ( const Point p,
const Error err 
)

Constructor for a fake vertex.

Definition at line 19 of file Vertex.cc.

References covariance_, dimension, i, idx(), j, and validity_.

19  :
20  chi2_( 0.0 ), ndof_( 0 ), position_( p ) {
21  index idx = 0;
22  for( index i = 0; i < dimension; ++ i )
23  for( index j = 0; j <= i; ++ j )
24  covariance_[ idx ++ ] = err( i, j );
25  validity_ = true;
26 }
int i
Definition: DBlmapReader.cc:9
float chi2_
chi-sqared
Definition: Vertex.h:153
unsigned int index
index type
Definition: Vertex.h:50
int j
Definition: DBlmapReader.cc:9
Point position_
position
Definition: Vertex.h:157
float covariance_[size]
covariance matrix (3x3) as vector
Definition: Vertex.h:159
index idx(index i, index j) const
position index
Definition: Vertex.h:170
float ndof_
number of degrees of freedom
Definition: Vertex.h:155
bool validity_
tells wether the vertex is really valid.
Definition: Vertex.h:166
Vertex::Vertex ( const Point p,
const Error err,
double  chi2,
double  ndof,
size_t  size 
)

constructor for a valid vertex, with all data

Definition at line 9 of file Vertex.cc.

References covariance_, dimension, i, idx(), j, tracks_, and validity_.

9  :
10  chi2_( chi2 ), ndof_( ndof ), position_( p ) {
11  tracks_.reserve( size );
12  index idx = 0;
13  for( index i = 0; i < dimension; ++ i )
14  for( index j = 0; j <= i; ++ j )
15  covariance_[ idx ++ ] = err( i, j );
16  validity_ = true;
17 }
int i
Definition: DBlmapReader.cc:9
std::vector< TrackBaseRef > tracks_
reference to tracks
Definition: Vertex.h:161
float chi2_
chi-sqared
Definition: Vertex.h:153
unsigned int index
index type
Definition: Vertex.h:50
double chi2() const
chi-squares
Definition: Vertex.h:82
int j
Definition: DBlmapReader.cc:9
double ndof() const
Definition: Vertex.h:89
Point position_
position
Definition: Vertex.h:157
float covariance_[size]
covariance matrix (3x3) as vector
Definition: Vertex.h:159
index idx(index i, index j) const
position index
Definition: Vertex.h:170
float ndof_
number of degrees of freedom
Definition: Vertex.h:155
bool validity_
tells wether the vertex is really valid.
Definition: Vertex.h:166

Member Function Documentation

void reco::Vertex::add ( const TrackBaseRef r,
float  w = 1.0 
)
void reco::Vertex::add ( const TrackBaseRef r,
const Track refTrack,
float  w = 1.0 
)

add the original a Track(reference) and the smoothed Track

double reco::Vertex::chi2 ( void  ) const
inline
double reco::Vertex::covariance ( int  i,
int  j 
) const
inline
CovarianceMatrix reco::Vertex::covariance ( void  ) const
inline

return SMatrix

Definition at line 114 of file Vertex.h.

References fill(), and m.

Referenced by xError(), yError(), and zError().

114 { Error m; fill( m ); return m; }
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
void fill(CovarianceMatrix &v) const
fill SMatrix
Definition: Vertex.cc:28
Error reco::Vertex::error ( ) const
inline

return SMatrix

Definition at line 116 of file Vertex.h.

References fill(), and m.

Referenced by argparse.ArgumentParser::_get_option_tuples(), python.rootplot.argparse.ArgumentParser::_get_option_tuples(), argparse.ArgumentParser::_parse_known_args(), python.rootplot.argparse.ArgumentParser::_parse_known_args(), argparse.ArgumentParser::_parse_optional(), python.rootplot.argparse.ArgumentParser::_parse_optional(), argparse.ArgumentParser::_read_args_from_files(), python.rootplot.argparse.ArgumentParser::_read_args_from_files(), IPTools::absoluteImpactParameter(), argparse.ArgumentParser::add_subparsers(), python.rootplot.argparse.ArgumentParser::add_subparsers(), FWSecVertexProxyBuilder::build(), FWVertexProxyBuilder::build(), reco::SecondaryVertex::computeDist2d(), reco::SecondaryVertex::computeDist3d(), dummyPrediction(), HLTDisplacedmumuFilter::hltFilter(), HLTDisplacedmumumuFilter::hltFilter(), pat::VertexAssociationSelector::operator()(), argparse.ArgumentParser::parse_args(), python.rootplot.argparse.ArgumentParser::parse_args(), argparse.ArgumentParser::parse_known_args(), python.rootplot.argparse.ArgumentParser::parse_known_args(), TrackIPProducer::produce(), ConeIsolation::produce(), SecondaryVertexProducer::produce(), and reco::GhostTrackVertexFinder::vertices().

116 { Error m; fill( m ); return m; }
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
void fill(CovarianceMatrix &v) const
fill SMatrix
Definition: Vertex.cc:28
void Vertex::fill ( CovarianceMatrix v) const

fill SMatrix

Definition at line 28 of file Vertex.cc.

References covariance_, dimension, i, idx(), and j.

Referenced by covariance(), error(), CandCommonVertexFitterBase::set(), and PFCandCommonVertexFitterBase::set().

28  {
29  index idx = 0;
30  for( index i = 0; i < dimension; ++ i )
31  for( index j = 0; j <= i; ++ j )
32  err( i, j ) = covariance_[ idx ++ ];
33 }
int i
Definition: DBlmapReader.cc:9
unsigned int index
index type
Definition: Vertex.h:50
int j
Definition: DBlmapReader.cc:9
float covariance_[size]
covariance matrix (3x3) as vector
Definition: Vertex.h:159
index idx(index i, index j) const
position index
Definition: Vertex.h:170
bool reco::Vertex::hasRefittedTracks ( ) const
inline

Checks whether refitted tracks are stored.

Definition at line 121 of file Vertex.h.

References refittedTracks_.

Referenced by nTracks(), reco::VertexFilter::operator()(), GhostTrackComputer::operator()(), CombinedSVComputer::operator()(), p4(), and reco::TrackKinematics::TrackKinematics().

121 {return !refittedTracks_.empty();}
std::vector< Track > refittedTracks_
The vector of refitted tracks.
Definition: Vertex.h:163
index reco::Vertex::idx ( index  i,
index  j 
) const
inlineprivate

position index

Definition at line 170 of file Vertex.h.

References a, b, and j.

Referenced by covariance(), fill(), and Vertex().

170  {
171  int a = ( i <= j ? i : j ), b = ( i <= j ? j : i );
172  return b * ( b + 1 ) / 2 + a;
173  }
int i
Definition: DBlmapReader.cc:9
int j
Definition: DBlmapReader.cc:9
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
bool reco::Vertex::isFake ( ) const
inline

Tells whether a Vertex is fake, i.e. not a vertex made out of a proper fit with tracks. For a primary vertex, it could simply be the beam line.

Definition at line 65 of file Vertex.h.

References chi2_, ndof_, and tracks_.

Referenced by TopDiLeptonDQM::analyze(), TopHLTDiMuonDQM::analyze(), PileupJetIdAlgo::computeIdVariables(), BVertexFilter::filter(), TauTagTools::filteredPFChargedHadrCands(), TauTagTools::filteredTracks(), HLTVertexFilter::hltFilter(), PVObjectSelector::operator()(), EventVtxInfoNtupleDumper::produce(), PFchsMETcorrInputProducer::produce(), and PrimaryVertexMonitor::vertexPlots().

65 {return (chi2_==0 && ndof_==0 && tracks_.empty());}
std::vector< TrackBaseRef > tracks_
reference to tracks
Definition: Vertex.h:161
float chi2_
chi-sqared
Definition: Vertex.h:153
float ndof_
number of degrees of freedom
Definition: Vertex.h:155
bool reco::Vertex::isValid ( void  ) const
inline
double reco::Vertex::ndof ( ) const
inline

Number of degrees of freedom Meant to be Double32_t for soft-assignment fitters: tracks may contribute to the vertex with fractional weights. The ndof is then = to the sum of the track weights. see e.g. CMS NOTE-2006/032, CMS NOTE-2004/002

Definition at line 89 of file Vertex.h.

References ndof_.

Referenced by CMSDAS11DijetAnalyzer::analyze(), CMSDAS11DijetTestAnalyzer::analyze(), VertexMonitor::analyze(), SiPixelTrackResidualSource::analyze(), TkConvValidator::analyze(), PrimaryVertexValidation::analyze(), ConversionProducer::buildCollection(), PileupJetIdAlgo::computeIdVariables(), reco::PFDisplacedVertex::Dump(), SVTagInfoValidationAnalyzer::fillRecoToSim(), recoBSVTagInfoValidationAnalyzer::fillRecoToSim(), recoBSVTagInfoValidationAnalyzer::fillSimToReco(), SVTagInfoValidationAnalyzer::fillSimToReco(), PrimaryVertexAnalyzer4PU::fillTrackHistos(), V0Fitter::fitAll(), PFDisplacedVertexFinder::fitVertexFromSeed(), HLTDisplacedmumuFilter::hltFilter(), HLTDisplacedmumumuFilter::hltFilter(), HLTVertexFilter::hltFilter(), ConversionTools::isGoodConversion(), PVObjectSelector::operator()(), EventVtxInfoNtupleDumper::produce(), ConeIsolation::produce(), PFchsMETcorrInputProducer::produce(), pat::PATConversionProducer::produce(), PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::selectPriVtxCompatibleWithTrack(), CandCommonVertexFitterBase::set(), PFCandCommonVertexFitterBase::set(), and PrimaryVertexMonitor::vertexPlots().

89 { return ndof_; }
float ndof_
number of degrees of freedom
Definition: Vertex.h:155
double reco::Vertex::normalizedChi2 ( ) const
inline
unsigned int Vertex::nTracks ( float  minWeight = 0.5) const

Returns the number of tracks in the vertex with weight above minWeight.

Definition at line 147 of file Vertex.cc.

References hasRefittedTracks(), n, originalTrack(), refittedTracks_, tracks_begin(), tracks_end(), and trackWeight().

148 {
149  int n=0;
150  if(hasRefittedTracks()) {
151  for(std::vector<Track>::const_iterator iter = refittedTracks_.begin(); iter != refittedTracks_.end(); ++iter)
152  if (trackWeight(originalTrack(*iter)) >=minWeight)
153  n++;
154  }
155  else
156  {
157  for(std::vector<reco::TrackBaseRef>::const_iterator iter = tracks_begin(); iter != tracks_end(); iter++)
158  if (trackWeight(*iter) >=minWeight)
159  n++;
160  }
161  return n;
162 }
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:45
TrackBaseRef originalTrack(const Track &refTrack) const
Definition: Vertex.cc:86
std::vector< Track > refittedTracks_
The vector of refitted tracks.
Definition: Vertex.h:163
bool hasRefittedTracks() const
Checks whether refitted tracks are stored.
Definition: Vertex.h:121
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:40
TrackBaseRef Vertex::originalTrack ( const Track 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 86 of file Vertex.cc.

References pos, refittedTracks_, and tracks_.

Referenced by PFDisplacedVertexFinder::commonTracks(), PF_PU_AssoMapAlgos::FindNIVertex(), PFDisplacedVertexHelper::identifyVertex(), PFDisplacedVertexHelper::lambdaCP(), reco::PFDisplacedVertex::momentum(), nTracks(), and p4().

87 {
88  if (refittedTracks_.empty())
89  throw cms::Exception("Vertex") << "No refitted tracks stored in vertex\n";
90  std::vector<Track>::const_iterator it =
91  find_if(refittedTracks_.begin(), refittedTracks_.end(), TrackEqual(refTrack));
92  if (it==refittedTracks_.end())
93  throw cms::Exception("Vertex") << "Refitted track not found in list\n";
94  size_t pos = it - refittedTracks_.begin();
95  return tracks_[pos];
96 }
std::vector< TrackBaseRef > tracks_
reference to tracks
Definition: Vertex.h:161
std::vector< Track > refittedTracks_
The vector of refitted tracks.
Definition: Vertex.h:163
math::XYZTLorentzVectorD Vertex::p4 ( float  mass = 0.13957018,
float  minWeight = 0.5 
) const

Returns the four momentum of the sum of the tracks, assuming the given mass for the decay products.

Definition at line 113 of file Vertex.cc.

References hasRefittedTracks(), originalTrack(), refittedTracks_, tracks_begin(), tracks_end(), and trackWeight().

Referenced by FWVertexProxyBuilder::build(), SecondaryVertexProducer::produce(), reco::Conversion::refittedPair4Momentum(), and TrackVertexArbitration::trackVertexArbitrator().

114 {
115 
117  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec;
118 
119  if(hasRefittedTracks()) {
120  for(std::vector<Track>::const_iterator iter = refittedTracks_.begin();
121  iter != refittedTracks_.end(); ++iter) {
122  if (trackWeight(originalTrack(*iter)) >=minWeight) {
123  vec.SetPx(iter->px());
124  vec.SetPy(iter->py());
125  vec.SetPz(iter->pz());
126  vec.SetM(mass);
127  sum += vec;
128  }
129  }
130  }
131  else
132  {
133  for(std::vector<reco::TrackBaseRef>::const_iterator iter = tracks_begin();
134  iter != tracks_end(); iter++) {
135  if (trackWeight(*iter) >=minWeight) {
136  vec.SetPx((*iter)->px());
137  vec.SetPy((*iter)->py());
138  vec.SetPz((*iter)->pz());
139  vec.SetM(mass);
140  sum += vec;
141  }
142  }
143  }
144  return sum;
145 }
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:15
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:45
TrackBaseRef originalTrack(const Track &refTrack) const
Definition: Vertex.cc:86
std::vector< Track > refittedTracks_
The vector of refitted tracks.
Definition: Vertex.h:163
bool hasRefittedTracks() const
Checks whether refitted tracks are stored.
Definition: Vertex.h:121
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:40
const Point& reco::Vertex::position ( ) const
inline

position

Definition at line 93 of file Vertex.h.

References position_.

Referenced by IPTools::absoluteImpactParameter(), IPTools::absoluteImpactParameter3D(), IPTools::absoluteTransverseImpactParameter(), CMSDAS11DijetAnalyzer::analyze(), CMSDAS11DijetTestAnalyzer::analyze(), TopElectronHLTOfflineSource::analyze(), TkConvValidator::analyze(), PrimaryVertexValidation::analyze(), SignedImpactParameter3D::apply(), CaloRecoTauAlgorithm::buildCaloTau(), PFRecoTauAlgorithm::buildPFTau(), pat::LeptonVertexSignificance::calculate(), ConversionProducer::checkPhi(), VertexDistance::compatibility(), reco::SecondaryVertex::computeDist2d(), reco::SecondaryVertex::computeDist3d(), PileupJetIdAlgo::computeIdVariables(), GsfConstraintAtVertex::constrainAtVertex(), reco::VZero::crossingPoint(), VertexDistance::distance(), VertexCompatibleWithBeam::distanceToBeam(), dummyPrediction(), reco::PFDisplacedVertex::Dump(), GsfVertexTrackCompatibilityEstimator::estimate(), KalmanVertexTrackCompatibilityEstimator< N >::estimate(), PrimaryVertexAnalyzer4PU::fillTrackHistos(), TrackingFailureFilter::filter(), TauTagTools::filteredPFChargedHadrCands(), TauTagTools::filteredTracks(), PF_PU_AssoMapAlgos::FindConversionVertex(), PF_PU_AssoMapAlgos::FindNIVertex(), PFDisplacedVertexFinder::fitVertexFromSeed(), BtoCharmDecayVertexMerger::flightDirection(), HLTDisplacedmumuFilter::hltFilter(), HLTDisplacedmumumuFilter::hltFilter(), PFDisplacedVertexHelper::identifyVertex(), NuclearVertexBuilder::isCompatible(), muon::isHighPtMuon(), PFTauQualityCutWrapper::isolationChargedObjects(), PFTauQualityCutWrapper::isolationPUObjects(), SoftPFMuonTagInfoProducer::isSoftMuon(), muon::isSoftMuon(), SoftPFMuonTagInfoProducer::isTightMuon(), muon::isTightMuon(), PhotonConversionTrajectorySeedProducerFromQuadrupletsAlgo::loop(), PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::loopOnPriVtx(), ConversionTools::matchesConversion(), TracksClusteringFromDisplacedSeed::nearTracks(), NuclearTrackCorrector::newTrajNeeded(), reco::VertexFilter::operator()(), VertexCompatibleWithBeam::operator()(), reco::tau::RecoTauEnergyRecoveryPlugin::operator()(), PVObjectSelector::operator()(), pat::VertexAssociationSelector::operator()(), print(), TrackIPProducer::produce(), EventVtxInfoNtupleDumper::produce(), TrackFromPVSelector::produce(), MuonFromPVSelector::produce(), InclusiveVertexFinder::produce(), GsfElectronFromPVSelector::produce(), VertexFromTrackProducer::produce(), SecondaryVertexProducer::produce(), TkConvValidator::recalculateMomentumAtFittedVertex(), pf2pat::IPCutPFCandidateSelectorDefinition::select(), PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::selectPriVtxCompatibleWithTrack(), reco::PFDisplacedVertex::setPrimaryDirection(), FWConvTrackHitsDetailView::setTextInfo(), reco::V0Candidate::setVertex(), PFTauQualityCutWrapper::signalChargedObjects(), VertexDistanceXY::signedDistance(), VertexDistance3D::signedDistance(), IPTools::signedImpactParameter3D(), IPTools::signedTransverseImpactParameter(), ConeIsolationAlgorithm::tag(), ConversionProducer::trackD0Cut(), QcdUeDQM::trackSelection(), TrackVertexArbitration::trackVertexArbitrator(), reco::JetSignalVertexCompatibilityAlgo::trackVertexCompat(), PrimaryVertexMonitor::vertexPlots(), and reco::GhostTrackVertexFinder::vertices().

93 { return position_; }
Point position_
position
Definition: Vertex.h:157
Track reco::Vertex::refittedTrack ( const TrackBaseRef 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

Referenced by CombinedSVComputer::operator()(), GhostTrackComputer::operator()(), refittedTrack(), BtoCharmDecayVertexMerger::resolveBtoDchain(), reco::TrackKinematics::TrackKinematics(), and reco::PFDisplacedVertex::trackPosition().

Track Vertex::refittedTrack ( const TrackRef 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 108 of file Vertex.cc.

References refittedTrack().

109 {
110  return refittedTrack(TrackBaseRef(track));
111 }
Track refittedTrack(const TrackBaseRef &track) const
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:22
const std::vector<Track>& reco::Vertex::refittedTracks ( ) const
inline
void Vertex::removeTracks ( )

Definition at line 65 of file Vertex.cc.

References refittedTracks_, tracks_, and weights_.

Referenced by reco::PFDisplacedVertex::cleanTracks(), and PFDisplacedVertexFinder::fitVertexFromSeed().

66 {
67  weights_.clear();
68  tracks_.clear();
69  refittedTracks_.clear();
70 }
std::vector< TrackBaseRef > tracks_
reference to tracks
Definition: Vertex.h:161
std::vector< Track > refittedTracks_
The vector of refitted tracks.
Definition: Vertex.h:163
std::vector< uint8_t > weights_
Definition: Vertex.h:164
Vertex::trackRef_iterator Vertex::tracks_begin ( ) const

first iterator over tracks

Definition at line 40 of file Vertex.cc.

References tracks_.

Referenced by NuclearVertexBuilder::addSecondaryTrack(), VertexMonitor::analyze(), PrimaryVertexValidation::analyze(), PrimaryVertexAnalyzer4PU::analyzeVertexCollectionTP(), VertexAssociatorByTracks::associateRecoToSim(), VertexAssociatorByTracks::associateSimToReco(), FWSecVertexProxyBuilder::build(), FWVertexProxyBuilder::build(), NuclearLikelihood::calculate(), PFPileUpAlgo::chargedHadronVertex(), PileupJetIdAlgo::computeIdVariables(), computeSharedTracks(), PFRecoTauDiscriminationByFlight::discriminate(), GsfVertexTrackCompatibilityEstimator::estimate(), KalmanVertexTrackCompatibilityEstimator< N >::estimate(), SVTagInfoValidationAnalyzer::fillRecoToSim(), recoBSVTagInfoValidationAnalyzer::fillRecoToSim(), SVTagInfoValidationAnalyzer::fillSimToReco(), recoBSVTagInfoValidationAnalyzer::fillSimToReco(), FFTJetPFPileupCleaner::findSomeVertexWFakes(), PrimaryVertexAnalyzer4PU::getTruthMatchedVertexTracks(), PATPrimaryVertexSelector::getVertexVariables(), HLTDisplacedmumumuFilter::hltFilter(), HLTDisplacedmumuFilter::hltFilter(), NuclearVertexBuilder::isCompatible(), main(), reco::VZero::negativeDaughter(), nTracks(), reco::VertexFilter::operator()(), p4(), reco::VZero::positiveDaughter(), reco::NuclearInteraction::primaryTrack(), PrimaryVertexAnalyzer4PU::printEventSummary(), PFchsMETcorrInputProducer::produce(), SecondaryVertexProducer::produce(), PVClusterComparer::pTSquaredSum(), BtoCharmDecayVertexMerger::resolveBtoDchain(), NuclearLikelihood::secondaryTrackMaxHits(), reco::NuclearInteraction::secondaryTracks_begin(), PFElectronAlgo::SetIDOutputs(), PFElectronAlgo::SetLinks(), PFEGammaAlgo::SetLinks(), reco::NbSharedTracks::sharedTracks(), VertexHigherPtSquared::sumPtSquared(), TrackerDpgAnalysis::sumPtSquared(), reco::TrackKinematics::TrackKinematics(), and PrimaryVertexMonitor::vertexPlots().

41 {
42  return tracks_.begin();
43 }
std::vector< TrackBaseRef > tracks_
reference to tracks
Definition: Vertex.h:161
Vertex::trackRef_iterator Vertex::tracks_end ( ) const

last iterator over tracks

Definition at line 45 of file Vertex.cc.

References tracks_.

Referenced by NuclearVertexBuilder::addSecondaryTrack(), VertexMonitor::analyze(), PrimaryVertexValidation::analyze(), PrimaryVertexAnalyzer4PU::analyzeVertexCollectionTP(), VertexAssociatorByTracks::associateRecoToSim(), VertexAssociatorByTracks::associateSimToReco(), FWSecVertexProxyBuilder::build(), FWVertexProxyBuilder::build(), PFPileUpAlgo::chargedHadronVertex(), PileupJetIdAlgo::computeIdVariables(), computeSharedTracks(), PFRecoTauDiscriminationByFlight::discriminate(), GsfVertexTrackCompatibilityEstimator::estimate(), KalmanVertexTrackCompatibilityEstimator< N >::estimate(), SVTagInfoValidationAnalyzer::fillRecoToSim(), recoBSVTagInfoValidationAnalyzer::fillRecoToSim(), SVTagInfoValidationAnalyzer::fillSimToReco(), recoBSVTagInfoValidationAnalyzer::fillSimToReco(), FFTJetPFPileupCleaner::findSomeVertexWFakes(), PrimaryVertexAnalyzer4PU::getTruthMatchedVertexTracks(), PATPrimaryVertexSelector::getVertexVariables(), main(), nTracks(), reco::VertexFilter::operator()(), p4(), PrimaryVertexAnalyzer4PU::printEventSummary(), PFchsMETcorrInputProducer::produce(), SecondaryVertexProducer::produce(), PVClusterComparer::pTSquaredSum(), BtoCharmDecayVertexMerger::resolveBtoDchain(), NuclearLikelihood::secondaryTrackMaxHits(), reco::NuclearInteraction::secondaryTracks_end(), PFElectronAlgo::SetIDOutputs(), PFElectronAlgo::SetLinks(), PFEGammaAlgo::SetLinks(), reco::NbSharedTracks::sharedTracks(), VertexHigherPtSquared::sumPtSquared(), TrackerDpgAnalysis::sumPtSquared(), reco::TrackKinematics::TrackKinematics(), and PrimaryVertexMonitor::vertexPlots().

46 {
47 // if ( !(tracks_.size() ) ) createTracks();
48  return tracks_.end();
49  // return weights_.keys().end();
50 }
std::vector< TrackBaseRef > tracks_
reference to tracks
Definition: Vertex.h:161
size_t Vertex::tracksSize ( ) const
float reco::Vertex::trackWeight ( const TrackBaseRef r) const
float Vertex::trackWeight ( const TrackRef r) const

returns the weight with which a Track has contributed to the vertex-fit.

Definition at line 80 of file Vertex.cc.

References trackWeight().

81 {
82  return trackWeight(TrackBaseRef(track));
83 }
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:22
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
double reco::Vertex::x ( ) const
inline

x coordinate

Definition at line 95 of file Vertex.h.

References position_.

Referenced by svgfig.Curve.Sample::__repr__(), svgfig.Ellipse::__repr__(), TopDiLeptonDQM::analyze(), TkConvValidator::analyze(), TopHLTDiMuonDQM::analyze(), SignedImpactParameter3D::apply(), SignedTransverseImpactParameter::apply(), SignedDecayLength3D::apply(), Vispa.Gui.WidgetContainer.WidgetContainer::autosize(), TauDiscriminationAgainstCaloMuon< TauType, TauDiscriminator >::beginEvent(), Vispa.Gui.VispaWidget.VispaWidget::boundingRect(), FWSecVertexProxyBuilder::build(), FWVertexProxyBuilder::build(), PFRecoTauAlgorithm::buildPFTau(), SignedDecayLength3D::closestApproachToJet(), SignedImpactParameter3D::closestApproachToJet(), IPTools::closestApproachToJet(), SignedImpactParameter3D::distance(), SignedImpactParameter3D::distanceWithJetAxis(), reco::Conversion::dxy(), reco::Conversion::dz(), PFPhotonAlgo::EvaluateSingleLegMVA(), PFEGammaAlgo::EvaluateSingleLegMVA(), SVTagInfoValidationAnalyzer::fillRecoToSim(), recoBSVTagInfoValidationAnalyzer::fillRecoToSim(), SVTagInfoValidationAnalyzer::fillSimToReco(), recoBSVTagInfoValidationAnalyzer::fillSimToReco(), V0Fitter::fitAll(), getVertexD0(), HLTVertexFilter::hltFilter(), ConversionTools::isGoodConversion(), IPTools::jetTrackDistance(), IPTools::linearizedSignedImpactParameter3D(), reco::Conversion::lxy(), reco::PFDisplacedVertex::momentum(), ConversionHitChecker::nHitsBeforeVtx(), pat::VertexAssociationSelector::operator()(), geometryXMLparser.Alignable::pos(), Vispa.Gui.ConnectableWidget.ConnectableWidget::positionizeMenuWidget(), reco::print(), ParticleReplacerParticleGun::produce(), pat::PATConversionProducer::produce(), SecondaryVertexProducer::produce(), PFAlgo::reconstructCluster(), PFEGammaAlgo::RunPFEG(), PFPhotonAlgo::RunPFPhoton(), reco::PFDisplacedVertex::setPrimaryDirection(), FWConvTrackHitsDetailView::setTextInfo(), IPTools::signedDecayLength3D(), IPTools::signedImpactParameter3D(), IPTools::signedTransverseImpactParameter(), TrackClusterSplitter::splitClusters(), and SignedTransverseImpactParameter::zImpactParameter().

95 { return position_.X(); }
Point position_
position
Definition: Vertex.h:157
double reco::Vertex::xError ( ) const
inline
double reco::Vertex::y ( ) const
inline

y coordinate

Definition at line 97 of file Vertex.h.

References position_.

Referenced by svgfig.Ellipse::__repr__(), TopDiLeptonDQM::analyze(), TkConvValidator::analyze(), TopHLTDiMuonDQM::analyze(), SignedImpactParameter3D::apply(), SignedTransverseImpactParameter::apply(), SignedDecayLength3D::apply(), Vispa.Gui.WidgetContainer.WidgetContainer::autosize(), TauDiscriminationAgainstCaloMuon< TauType, TauDiscriminator >::beginEvent(), Vispa.Gui.VispaWidget.VispaWidget::boundingRect(), FWSecVertexProxyBuilder::build(), FWVertexProxyBuilder::build(), PFRecoTauAlgorithm::buildPFTau(), SignedDecayLength3D::closestApproachToJet(), SignedImpactParameter3D::closestApproachToJet(), IPTools::closestApproachToJet(), SignedImpactParameter3D::distance(), SignedImpactParameter3D::distanceWithJetAxis(), reco::Conversion::dxy(), reco::Conversion::dz(), PFPhotonAlgo::EvaluateSingleLegMVA(), PFEGammaAlgo::EvaluateSingleLegMVA(), recoBSVTagInfoValidationAnalyzer::fillRecoToSim(), SVTagInfoValidationAnalyzer::fillRecoToSim(), SVTagInfoValidationAnalyzer::fillSimToReco(), recoBSVTagInfoValidationAnalyzer::fillSimToReco(), V0Fitter::fitAll(), getVertexD0(), HLTVertexFilter::hltFilter(), ConversionTools::isGoodConversion(), IPTools::jetTrackDistance(), IPTools::linearizedSignedImpactParameter3D(), reco::Conversion::lxy(), reco::PFDisplacedVertex::momentum(), ConversionHitChecker::nHitsBeforeVtx(), pat::VertexAssociationSelector::operator()(), geometryXMLparser.Alignable::pos(), Vispa.Gui.ConnectableWidget.ConnectableWidget::positionizeMenuWidget(), reco::print(), ParticleReplacerParticleGun::produce(), pat::PATConversionProducer::produce(), SecondaryVertexProducer::produce(), PFAlgo::reconstructCluster(), PFEGammaAlgo::RunPFEG(), PFPhotonAlgo::RunPFPhoton(), reco::PFDisplacedVertex::setPrimaryDirection(), FWConvTrackHitsDetailView::setTextInfo(), IPTools::signedDecayLength3D(), IPTools::signedImpactParameter3D(), IPTools::signedTransverseImpactParameter(), TrackClusterSplitter::splitClusters(), and SignedTransverseImpactParameter::zImpactParameter().

97 { return position_.Y(); }
Point position_
position
Definition: Vertex.h:157
double reco::Vertex::yError ( ) const
inline
double reco::Vertex::z ( ) const
inline

y coordinate

Definition at line 99 of file Vertex.h.

References position_.

Referenced by CMSDAS11DijetAnalyzer::analyze(), CMSDAS11DijetTestAnalyzer::analyze(), TopDiLeptonDQM::analyze(), TopHLTDiMuonDQM::analyze(), SignedImpactParameter3D::apply(), SignedTransverseImpactParameter::apply(), SignedDecayLength3D::apply(), TauDiscriminationAgainstCaloMuon< TauType, TauDiscriminator >::beginEvent(), FWSecVertexProxyBuilder::build(), FWVertexProxyBuilder::build(), CaloRecoTauTagInfoAlgorithm::buildCaloTauTagInfo(), PFRecoTauAlgorithm::buildPFTau(), PFRecoTauTagInfoAlgorithm::buildPFTauTagInfo(), SignedDecayLength3D::closestApproachToJet(), SignedImpactParameter3D::closestApproachToJet(), IPTools::closestApproachToJet(), SignedImpactParameter3D::distance(), SignedImpactParameter3D::distanceWithJetAxis(), reco::Conversion::dz(), PFPhotonAlgo::EvaluateLCorrMVA(), PFEGammaAlgo::EvaluateLCorrMVA(), PFPhotonAlgo::EvaluateSingleLegMVA(), PFEGammaAlgo::EvaluateSingleLegMVA(), EwkElecTauHistManager::fillHistograms(), EwkMuTauHistManager::fillHistograms(), recoBSVTagInfoValidationAnalyzer::fillRecoToSim(), SVTagInfoValidationAnalyzer::fillRecoToSim(), SVTagInfoValidationAnalyzer::fillSimToReco(), recoBSVTagInfoValidationAnalyzer::fillSimToReco(), V0Fitter::fitAll(), HLTVertexFilter::hltFilter(), IPTools::jetTrackDistance(), IPTools::linearizedSignedImpactParameter3D(), QualityCutsAnalyzer::LoopOverJetTracksAssociation(), reco::Conversion::lz(), PrimaryVertexAnalyzer::matchVertex(), PrimaryVertexAnalyzer4PU::matchVertex(), reco::PFDisplacedVertex::momentum(), ConversionHitChecker::nHitsBeforeVtx(), PVObjectSelector::operator()(), pat::VertexAssociationSelector::operator()(), geometryXMLparser.Alignable::pos(), reco::print(), EventVtxInfoNtupleDumper::produce(), ConeIsolation::produce(), ParticleReplacerParticleGun::produce(), PFchsMETcorrInputProducer::produce(), SecondaryVertexProducer::produce(), PFAlgo::reconstructCluster(), PFEGammaAlgo::RunPFEG(), PFPhotonAlgo::RunPFPhoton(), reco::PFDisplacedVertex::setPrimaryDirection(), FWConvTrackHitsDetailView::setTextInfo(), IPTools::signedDecayLength3D(), IPTools::signedImpactParameter3D(), TrackClusterSplitter::splitClusters(), ConeIsolationAlgorithm::tag(), and SignedTransverseImpactParameter::zImpactParameter().

99 { return position_.Z(); }
Point position_
position
Definition: Vertex.h:157
double reco::Vertex::zError ( ) const
inline

Member Data Documentation

float reco::Vertex::chi2_
private

chi-sqared

Definition at line 153 of file Vertex.h.

Referenced by chi2(), isFake(), and normalizedChi2().

float reco::Vertex::covariance_[size]
private

covariance matrix (3x3) as vector

Definition at line 159 of file Vertex.h.

Referenced by covariance(), fill(), and Vertex().

float reco::Vertex::ndof_
private

number of degrees of freedom

Definition at line 155 of file Vertex.h.

Referenced by isFake(), ndof(), and normalizedChi2().

Point reco::Vertex::position_
private

position

Definition at line 157 of file Vertex.h.

Referenced by position(), x(), y(), and z().

std::vector<Track> reco::Vertex::refittedTracks_
private

The vector of refitted tracks.

Definition at line 163 of file Vertex.h.

Referenced by hasRefittedTracks(), nTracks(), originalTrack(), p4(), refittedTracks(), and removeTracks().

std::vector<TrackBaseRef > reco::Vertex::tracks_
private

reference to tracks

Definition at line 161 of file Vertex.h.

Referenced by isFake(), originalTrack(), removeTracks(), tracks_begin(), tracks_end(), and Vertex().

bool reco::Vertex::validity_
private

tells wether the vertex is really valid.

Definition at line 166 of file Vertex.h.

Referenced by isValid(), and Vertex().

std::vector<uint8_t> reco::Vertex::weights_
private

Definition at line 164 of file Vertex.h.

Referenced by removeTracks(), and tracksSize().