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 | Static Public Member Functions | Static Public Attributes | Private Attributes
reco::TrackBase Class Reference

#include <DataFormats/TrackReco/interface/TrackBase.h>

Inheritance diagram for reco::TrackBase:
reco::Track reco::GsfTrack reco::TrackTransientTrack reco::GsfTransientTrack

Public Types

enum  { dimension = 5 }
 parameter dimension More...
 
enum  { covarianceSize = dimension * ( dimension + 1 ) / 2 }
 error matrix size More...
 
enum  {
  i_qoverp = 0, i_lambda, i_phi, i_dxy,
  i_dsz
}
 enumerator provided indices to the five parameters More...
 
typedef math::Error< dimension >
::type 
CovarianceMatrix
 5 parameter covariance matrix More...
 
typedef unsigned int index
 index type More...
 
typedef math::Vector
< dimension >::type 
ParameterVector
 parameter vector More...
 
typedef math::XYZPoint Point
 point in the space More...
 
enum  TrackAlgorithm {
  undefAlgorithm =0, ctf =1, rs =2, cosmics =3,
  iter0 =4, iter1 =5, iter2 =6, iter3 =7,
  iter4 =8, iter5 =9, iter6 =10, iter7 =11,
  iter8 =12, iter9 =13, iter10 =14, outInEcalSeededConv =15,
  inOutEcalSeededConv =16, nuclInter =17, standAloneMuon =18, globalMuon =19,
  cosmicStandAloneMuon =20, cosmicGlobalMuon =21, iter1LargeD0 =22, iter2LargeD0 =23,
  iter3LargeD0 =24, iter4LargeD0 =25, iter5LargeD0 =26, bTagGhostTracks =27,
  beamhalo =28, gsf =29, algoSize =30
}
 track algorithm More...
 
enum  TrackQuality {
  undefQuality =-1, loose =0, tight =1, highPurity =2,
  confirmed =3, goodIterative =4, looseSetWithPV =5, highPuritySetWithPV =6,
  qualitySize =7
}
 track quality More...
 
typedef math::XYZVector Vector
 spatial vector More...
 

Public Member Functions

TrackAlgorithm algo () const
 
std::string algoName () const
 
void appendHitPattern (const TrackingRecHit &hit)
 
int charge () const
 track electric charge More...
 
double chi2 () const
 chi-squared of the fit More...
 
CovarianceMatrix covariance () const
 return track covariance matrix More...
 
double covariance (int i, int j) const
 (i,j)-th element of covarianve matrix ( i, j = 0, ... 4 ) More...
 
double d0 () const
 dxy parameter in perigee convention (d0 = - dxy) More...
 
double d0Error () const
 error on d0 More...
 
double dsz () const
 dsz parameter (THIS IS NOT the SZ impact parameter to (0,0,0) if refPoint is far from (0,0,0): see parametrization definition above for details) More...
 
double dsz (const Point &myBeamSpot) const
 dsz parameter with respect to a user-given beamSpot (WARNING: this quantity can only be interpreted as the distance in the S-Z plane to the beamSpot, if the beam spot is reasonably close to the refPoint, since linear approximations are involved). This is a good approximation for Tracker tracks. More...
 
double dszError () const
 error on dsz More...
 
double dxy () const
 dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close to (0,0,0): see parametrization definition above for details). See also function dxy(myBeamSpot) below. More...
 
double dxy (const Point &myBeamSpot) const
 dxy parameter with respect to a user-given beamSpot (WARNING: this quantity can only be interpreted as a minimum transverse distance if beamSpot, if the beam spot is reasonably close to the refPoint, since linear approximations are involved). This is a good approximation for Tracker tracks. More...
 
double dxy (const BeamSpot &theBeamSpot) const
 dxy parameter with respect to the beamSpot taking into account the beamspot slopes WARNING: this quantity can only be interpreted as a minimum transverse distance if beamSpot, if the beam spot is reasonably close to the refPoint, since linear approximations are involved). This is a good approximation for Tracker tracks. More...
 
double dxyError () const
 error on dxy More...
 
double dz () const
 dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to (0,0,0). See also function dz(myBeamSpot) below. More...
 
double dz (const Point &myBeamSpot) const
 dz parameter with respect to a user-given beamSpot (WARNING: this quantity can only be interpreted as the track z0, if the beamSpot is reasonably close to the refPoint, since linear approximations are involved). This is a good approximation for Tracker tracks. More...
 
double dzError () const
 error on dz More...
 
double error (int i) const
 error on specified element More...
 
double eta () const
 pseudorapidity of momentum vector More...
 
double etaError () const
 error on eta More...
 
CovarianceMatrixfill (CovarianceMatrix &v) const
 fill SMatrix More...
 
const HitPatternhitPattern () const
 Access the hit pattern, indicating in which Tracker layers the track has hits. More...
 
bool isLooper () const
 
double lambda () const
 Lambda angle. More...
 
double lambdaError () const
 error on lambda More...
 
const Vectormomentum () const
 track momentum vector More...
 
double ndof () const
 number of degrees of freedom of the fit More...
 
signed char nLoops () const
 
double normalizedChi2 () const
 chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero) More...
 
unsigned short numberOfLostHits () const
 number of cases where track crossed a layer without getting a hit. More...
 
unsigned short numberOfValidHits () const
 number of valid hits found More...
 
double p () const
 momentum vector magnitude More...
 
double parameter (int i) const
 i-th parameter ( i = 0, ... 4 ) More...
 
ParameterVector parameters () const
 Track parameters with one-to-one correspondence to the covariance matrix. More...
 
double phi () const
 azimuthal angle of momentum vector More...
 
double phiError () const
 error on phi More...
 
double pt () const
 track transverse momentum More...
 
double ptError () const
 error on Pt (set to 1000 TeV if charge==0 for safety) More...
 
double px () const
 x coordinate of momentum vector More...
 
double py () const
 y coordinate of momentum vector More...
 
double pz () const
 z coordinate of momentum vector More...
 
double qoverp () const
 q/p More...
 
double qoverpError () const
 error on signed transverse curvature More...
 
bool quality (const TrackQuality) const
 Track quality. More...
 
int qualityMask () const
 
const PointreferencePoint () const
 Reference point on the track. More...
 
void setAlgorithm (const TrackAlgorithm a, bool set=true)
 position index More...
 
template<typename C >
void setHitPattern (const C &c)
 set hit patterns from vector of hit references More...
 
template<typename I >
void setHitPattern (const I &begin, const I &end)
 
void setHitPattern (const TrackingRecHit &hit, size_t i)
 set hit pattern for specified hit More...
 
void setHitPattern (const HitPattern &hitP)
 set hitPattern from pre-defined hitPattern More...
 
void setNLoops (signed char value)
 
void setQuality (const TrackQuality, bool set=true)
 
void setQualityMask (int qualMask)
 
template<typename C >
void setTrackerExpectedHitsInner (const C &c)
 
template<typename I >
void setTrackerExpectedHitsInner (const I &begin, const I &end)
 
void setTrackerExpectedHitsInner (const TrackingRecHit &hit, size_t i)
 
void setTrackerExpectedHitsInner (const HitPattern &hitP)
 
template<typename C >
void setTrackerExpectedHitsOuter (const C &c)
 
template<typename I >
void setTrackerExpectedHitsOuter (const I &begin, const I &end)
 
void setTrackerExpectedHitsOuter (const TrackingRecHit &hit, size_t i)
 
void setTrackerExpectedHitsOuter (const HitPattern &hitP)
 
double theta () const
 polar angle More...
 
double thetaError () const
 error on theta More...
 
 TrackBase ()
 default constructor More...
 
 TrackBase (double chi2, double ndof, const Point &referencePoint, const Vector &momentum, int charge, const CovarianceMatrix &, TrackAlgorithm=undefAlgorithm, TrackQuality quality=undefQuality, signed char nloops=0)
 constructor from fit parameters and error matrix More...
 
const HitPatterntrackerExpectedHitsInner () const
 Access the hit pattern counting (in the Tracker) the number of expected crossed layers before the first trajectory's hit. More...
 
const HitPatterntrackerExpectedHitsOuter () const
 Access the hit pattern counting (in the Tracker) the number of expected crossed layers after the last trajectory's hit. More...
 
double validFraction () const
 fraction of valid hits on the track More...
 
const Pointvertex () const
 reference point on the track. This method is DEPRECATED, please use referencePoint() instead More...
 
double vx () const
 x coordinate of the reference point on track More...
 
double vy () const
 y coordinate of the reference point on track More...
 
double vz () const
 z coordinate of the reference point on track More...
 
 ~TrackBase ()
 virtual destructor More...
 

Static Public Member Functions

static TrackAlgorithm algoByName (const std::string &name)
 
static std::string algoName (TrackAlgorithm)
 
static index covIndex (index i, index j)
 covariance matrix index in array More...
 
static TrackQuality qualityByName (const std::string &name)
 
static std::string qualityName (TrackQuality)
 

Static Public Attributes

static const std::string algoNames []
 
static const std::string qualityNames [] = { "loose", "tight", "highPurity", "confirmed", "goodIterative", "looseSetWithPV", "highPuritySetWithPV"}
 

Private Attributes

uint8_t algorithm_
 track algorithm More...
 
char charge_
 electric charge More...
 
float chi2_
 chi-squared More...
 
float covariance_ [covarianceSize]
 perigee 5x5 covariance matrix More...
 
HitPattern hitPattern_
 hit pattern More...
 
Vector momentum_
 momentum vector at innermost point More...
 
float ndof_
 number of degrees of freedom More...
 
signed char nLoops_
 number of loops made during the building of the trajectory of a looper particle More...
 
uint8_t quality_
 track quality More...
 
HitPattern trackerExpectedHitsInner_
 hit pattern used for expected crossed layers after the last trajectory's hit More...
 
HitPattern trackerExpectedHitsOuter_
 hit pattern used for expected crossed layers before the first trajectory's hit More...
 
Point vertex_
 innermost (reference) point on track More...
 

Detailed Description

Common base class to all track types, including Muon fits. Internally, the following information is stored:

A reference position on the track: (vx,vy,vz)

Momentum at this given reference point on track: (px,py,pz)

5D curvilinear covariance matrix from the track fit

Charge

Chi-square and number of degrees of freedom

Summary information of the hit pattern

For tracks reconstructed in the CMS Tracker, the reference position is the point of closest approach to the centre of CMS. For muons, this is not necessarily true.

Parameters associated to the 5D curvilinear covariance matrix:
(qoverp, lambda, phi, dxy, dsz)
defined as:

qoverp = q / abs(p) = signed inverse of momentum [1/GeV]

lambda = pi/2 - polar angle at the given point

phi = azimuth angle at the given point

dxy = -vx*sin(phi) + vy*cos(phi) [cm]

dsz = vz*cos(lambda) - (vx*cos(phi)+vy*sin(phi))*sin(lambda) [cm]

Geometrically, dxy is the signed distance in the XY plane between the the straight line passing through (vx,vy) with azimuthal angle phi and the point (0,0).
The dsz parameter is the signed distance in the SZ plane between the the straight line passing through (vx,vy,vz) with angles (phi, lambda) and the point (s=0,z=0). The S axis is defined by the projection of the straight line onto the XY plane. The convention is to assign the S coordinate for (vx,vy) as the value vx*cos(phi)+vy*sin(phi). This value is zero when (vx,vy) is the point of minimum transverse distance to (0,0).

Note that dxy and dsz provide sensible estimates of the distance from the true particle trajectory to (0,0,0) ONLY in two cases:

When (vx,vy,vz) already correspond to the point of minimum transverse distance to (0,0,0) or is close to it (so that the differences between considering the exact trajectory or a straight line in this range are negligible). This is usually true for Tracker tracks.

When the track has infinite or extremely high momentum

More details about this parametrization are provided in the following document:
A. Strandlie, W. Wittek, "Propagation of Covariance Matrices...", CMS Note 2006/001

Author
Thomas Speer, Luca Lista, Pascal Vanlaer, Juan Alcaraz
Version
Id:
TrackBase.h,v 1.87 2012/04/26 10:12:39 bellan Exp

Definition at line 63 of file TrackBase.h.

Member Typedef Documentation

5 parameter covariance matrix

Definition at line 72 of file TrackBase.h.

typedef unsigned int reco::TrackBase::index

index type

Definition at line 80 of file TrackBase.h.

parameter vector

Definition at line 70 of file TrackBase.h.

point in the space

Definition at line 76 of file TrackBase.h.

spatial vector

Definition at line 74 of file TrackBase.h.

Member Enumeration Documentation

anonymous enum

parameter dimension

Enumerator
dimension 

Definition at line 66 of file TrackBase.h.

anonymous enum

error matrix size

Enumerator
covarianceSize 

Definition at line 68 of file TrackBase.h.

anonymous enum

enumerator provided indices to the five parameters

Enumerator
i_qoverp 
i_lambda 
i_phi 
i_dxy 
i_dsz 

Definition at line 78 of file TrackBase.h.

track algorithm

Enumerator
undefAlgorithm 
ctf 
rs 
cosmics 
iter0 
iter1 
iter2 
iter3 
iter4 
iter5 
iter6 
iter7 
iter8 
iter9 
iter10 
outInEcalSeededConv 
inOutEcalSeededConv 
nuclInter 
standAloneMuon 
globalMuon 
cosmicStandAloneMuon 
cosmicGlobalMuon 
iter1LargeD0 
iter2LargeD0 
iter3LargeD0 
iter4LargeD0 
iter5LargeD0 
bTagGhostTracks 
beamhalo 
gsf 
algoSize 

Definition at line 82 of file TrackBase.h.

82  { undefAlgorithm=0, ctf=1, rs=2, cosmics=3, iter0=4,
83  iter1=5, iter2=6, iter3=7, iter4=8, iter5=9, iter6=10, iter7=11, iter8=12, iter9=13,iter10=14,
85  nuclInter=17,
88  bTagGhostTracks=27,
89  beamhalo=28,
90  gsf=29,
91  algoSize=30 };

track quality

Enumerator
undefQuality 
loose 
tight 
highPurity 
confirmed 
goodIterative 
looseSetWithPV 
highPuritySetWithPV 
qualitySize 

Definition at line 95 of file TrackBase.h.

Constructor & Destructor Documentation

TrackBase::TrackBase ( )

default constructor

Definition at line 20 of file TrackBase.cc.

References covariance_, dimension, i, and j.

20  :
21  chi2_(0), vertex_(0,0,0), momentum_(0,0,0), ndof_(0), charge_(0), algorithm_(undefAlgorithm), quality_(0), nLoops_(0) {
22  index idx = 0;
23  for( index i = 0; i < dimension; ++ i )
24  for( index j = 0; j <= i; ++ j )
25  covariance_[ idx ++ ]=0;
26 }
float chi2_
chi-squared
Definition: TrackBase.h:306
int i
Definition: DBlmapReader.cc:9
unsigned int index
index type
Definition: TrackBase.h:80
uint8_t quality_
track quality
Definition: TrackBase.h:320
char charge_
electric charge
Definition: TrackBase.h:316
Point vertex_
innermost (reference) point on track
Definition: TrackBase.h:308
signed char nLoops_
number of loops made during the building of the trajectory of a looper particle
Definition: TrackBase.h:323
int j
Definition: DBlmapReader.cc:9
float ndof_
number of degrees of freedom
Definition: TrackBase.h:313
float covariance_[covarianceSize]
perigee 5x5 covariance matrix
Definition: TrackBase.h:304
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:310
uint8_t algorithm_
track algorithm
Definition: TrackBase.h:318
TrackBase::TrackBase ( double  chi2,
double  ndof,
const Point referencePoint,
const Vector momentum,
int  charge,
const CovarianceMatrix cov,
TrackAlgorithm  algorithm = undefAlgorithm,
TrackQuality  quality = undefQuality,
signed char  nloops = 0 
)

constructor from fit parameters and error matrix

Definition at line 28 of file TrackBase.cc.

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

30  :
31  chi2_( chi2 ), vertex_( vertex ), momentum_( momentum ), ndof_( ndof ), charge_( charge ), algorithm_(algorithm), quality_(0), nLoops_(nloops) {
32  index idx = 0;
33  for( index i = 0; i < dimension; ++ i )
34  for( index j = 0; j <= i; ++ j )
35  covariance_[ idx ++ ] = cov( i, j );
36  setQuality(quality);
37 }
float chi2_
chi-squared
Definition: TrackBase.h:306
int i
Definition: DBlmapReader.cc:9
< trclass="colgroup">< tdclass="colgroup"colspan=5 > Ecal cluster collections</td ></tr >< tr >< td >< ahref="classreco_1_1BasicCluster.html"> reco::BasicCluster</a ></td >< td >< ahref="DataFormats_EgammaReco.html"> reco::BasicClusterCollection</a ></td >< td >< ahref="#"> hybridSuperClusters</a ></td >< tdclass="description"> Basic clusters reconstructed with hybrid algorithm(barrel only)</td >< td >S.Rahatlou</td ></tr >< tr >< td >< a href
unsigned int index
index type
Definition: TrackBase.h:80
uint8_t quality_
track quality
Definition: TrackBase.h:320
char charge_
electric charge
Definition: TrackBase.h:316
Point vertex_
innermost (reference) point on track
Definition: TrackBase.h:308
signed char nLoops_
number of loops made during the building of the trajectory of a looper particle
Definition: TrackBase.h:323
void setQuality(const TrackQuality, bool set=true)
Definition: TrackBase.h:393
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:107
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:109
int j
Definition: DBlmapReader.cc:9
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:156
float ndof_
number of degrees of freedom
Definition: TrackBase.h:313
float covariance_[covarianceSize]
perigee 5x5 covariance matrix
Definition: TrackBase.h:304
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:310
int charge() const
track electric charge
Definition: TrackBase.h:113
uint8_t algorithm_
track algorithm
Definition: TrackBase.h:318
TrackBase::~TrackBase ( )

virtual destructor

Definition at line 39 of file TrackBase.cc.

39  {
40 }

Member Function Documentation

TrackBase::TrackAlgorithm reco::TrackBase::algo ( ) const
inline
TrackBase::TrackAlgorithm TrackBase::algoByName ( const std::string &  name)
static

Definition at line 55 of file TrackBase.cc.

References algoNames, algoSize, spr::find(), findQualityFiles::size, and undefAlgorithm.

Referenced by AlignmentTrackSelector::AlignmentTrackSelector(), BeamFitter::BeamFitter(), QcdUeDQM::QcdUeDQM(), and RecoTrackSelector::RecoTrackSelector().

55  {
58  if(index == size) return undefAlgorithm; // better this or throw() ?
59 
60  // cast
61  return TrackAlgorithm(index);
62 }
unsigned int index
index type
Definition: TrackBase.h:80
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
TrackAlgorithm
track algorithm
Definition: TrackBase.h:82
static const std::string algoNames[]
Definition: TrackBase.h:92
tuple size
Write out results.
std::string reco::TrackBase::algoName ( ) const
inline

Definition at line 336 of file TrackBase.h.

References algorithm_, beamhalo, bTagGhostTracks, cosmicGlobalMuon, cosmics, cosmicStandAloneMuon, ctf, globalMuon, gsf, inOutEcalSeededConv, iter0, iter1, iter10, iter1LargeD0, iter2, iter2LargeD0, iter3, iter3LargeD0, iter4, iter4LargeD0, iter5, iter5LargeD0, iter6, iter7, iter8, iter9, nuclInter, outInEcalSeededConv, rs, standAloneMuon, and undefAlgorithm.

336  {
337  // I'd like to do:
338  // return TrackBase::algoName(algorithm_);
339  // but I cannot define a const static function. Why???
340 
341  switch(algorithm_)
342  {
343  case undefAlgorithm: return "undefAlgorithm";
344  case ctf: return "ctf";
345  case rs: return "rs";
346  case cosmics: return "cosmics";
347  case beamhalo: return "beamhalo";
348  case iter0: return "iter0";
349  case iter1: return "iter1";
350  case iter2: return "iter2";
351  case iter3: return "iter3";
352  case iter4: return "iter4";
353  case iter5: return "iter5";
354  case iter6: return "iter6";
355  case iter7: return "iter7";
356  case iter8: return "iter8";
357  case iter9: return "iter9";
358  case iter10: return "iter10";
359  case outInEcalSeededConv: return "outInEcalSeededConv";
360  case inOutEcalSeededConv: return "inOutEcalSeededConv";
361  case nuclInter: return "nuclInter";
362  case standAloneMuon: return "standAloneMuon";
363  case globalMuon: return "globalMuon";
364  case cosmicStandAloneMuon: return "cosmicStandAloneMuon";
365  case cosmicGlobalMuon: return "cosmicGlobalMuon";
366  case iter1LargeD0: return "iter1LargeD0";
367  case iter2LargeD0: return "iter2LargeD0";
368  case iter3LargeD0: return "iter3LargeD0";
369  case iter4LargeD0: return "iter4LargeD0";
370  case iter5LargeD0: return "iter5LargeD0";
371  case bTagGhostTracks: return "bTagGhostTracks";
372  case gsf: return "gsf";
373  }
374  return "undefAlgorithm";
375  }
uint8_t algorithm_
track algorithm
Definition: TrackBase.h:318
std::string reco::TrackBase::algoName ( TrackAlgorithm  a)
inlinestatic

Definition at line 411 of file TrackBase.h.

References algoNames, and algoSize.

411  {
412  if(int(a) < int(algoSize) && int(a)>0) return algoNames[int(a)];
413  return "undefAlgorithm";
414  }
static const std::string algoNames[]
Definition: TrackBase.h:92
double a
Definition: hdecay.h:121
void reco::TrackBase::appendHitPattern ( const TrackingRecHit hit)
inline

Definition at line 261 of file TrackBase.h.

References reco::HitPattern::appendHit(), and hitPattern_.

261 { hitPattern_.appendHit( hit); }
void appendHit(const TrackingRecHit &hit)
Definition: HitPattern.cc:88
HitPattern hitPattern_
hit pattern
Definition: TrackBase.h:297
int reco::TrackBase::charge ( ) const
inline

track electric charge

Definition at line 113 of file TrackBase.h.

References charge_.

Referenced by PFTrackTransformer::addPoints(), PFTrackTransformer::addPointsAndBrems(), MultiTrackValidator::analyze(), EwkMuLumiMonitorDQM::analyze(), ValidationMisalignedTracker::analyze(), TestOutliers::analyze(), CosmicSplitterValidation::analyze(), MuonTrackValidator::analyze(), PrimaryVertexAnalyzer4PU::analyzeVertexCollection(), AlignmentTrackSelector::basicCuts(), FWSecVertexProxyBuilder::build(), FWVertexProxyBuilder::build(), ConversionProducer::buildCollection(), FWTrackProxyBuilderFF::buildTrack(), reco::TrackTransientTrack::charge(), reco::GsfTransientTrack::charge(), Tau3MuReco::check4MuonTrack(), AlignmentTwoBodyDecayTrackSelector::checkCharge(), converter::StandAloneMuonTrackToCandidate::convert(), converter::TrackToCandidate::convert(), muonisolation::PixelTrackExtractor::directionAtPresetRadius(), IdealHelixParameters::evalCircleCenter(), MTVHistoProducerAlgoForTracker::fill_simAssociated_recoTrack_histos(), TrackAnalyzer::fillHistosForState(), TrackerValidationVariables::fillTrackQuantities(), CSCEfficiency::filter(), cms::HICFTSfromL1orL2::FTSfromStandAlone(), ConversionFinder::getConversionInfo(), TrackDetectorAssociator::getFreeTrajectoryState(), HTrackAssociator::getFreeTrajectoryState(), AlCaHOCalibProducer::getFreeTrajectoryState(), VZeroFinder::getGlobalTrajectoryParameters(), PFDisplacedVertexCandidateFinder::getGlobalTrajectoryParameters(), FWPFTrackUtils::getTrack(), ValidHitPairFilter::getTrajectory(), EcalShowerProperties::getTrajectoryAtOuterPoint(), PlotRecTracks::getTrajectoryAtOuterPoint(), SeedFromProtoTrack::init(), trajectoryStateTransform::initialFreeState(), HLTmmkFilter::initialFreeState(), HLTmmkkFilter::initialFreeState(), trajectoryStateTransform::innerFreeState(), PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::inspectTrack(), MuonIdProducer::makeMuon(), MuonErrorMatrixAdjuster::makeTrack(), MatcherByPullsAlgorithm::match(), muonisolation::CaloExtractor::MuonAtCaloPosition(), reco::tau::RecoTauPiZeroStripPlugin2::operator()(), trajectoryStateTransform::outerFreeState(), fireworks::prepareTrack(), PrintRecoObjects::print(), IsolatedTracksCone::printTrack(), IsolatedTracksNxN::printTrack(), AlignmentMonitorMuonVsCurvature::processMuonResidualsFromTrack(), AlignmentMonitorSegmentDifferences::processMuonResidualsFromTrack(), AlignmentMonitorMuonSystemMap1D::processMuonResidualsFromTrack(), MuonAlignmentFromReference::processMuonResidualsFromTrack(), ShallowTracksProducer::produce(), QuarkoniaTrackSelector::produce(), PFDisplacedTrackerVertexProducer::produce(), SiStripElectronAssociator::produce(), ConversionSeedFilter::produce(), spr::propagateECAL(), spr::propagateHCAL(), spr::propagateTracker(), spr::propagateTrackerEnd(), spr::propagateTrackToECAL(), spr::propagateTrackToHCAL(), TrackExtrapolator::propagateTrackToVolume(), ptError(), qoverp(), TrackClassifier::reconstructionInformation(), PFAlgo::reconstructTrack(), MuonDTLocalMillepedeAlgorithm::run(), FWTrackHitsDetailView::setTextInfo(), FWConvTrackHitsDetailView::setTextInfo(), and ConversionProducer::trackD0Cut().

113 { return charge_; }
char charge_
electric charge
Definition: TrackBase.h:316
double reco::TrackBase::chi2 ( void  ) const
inline
CovarianceMatrix reco::TrackBase::covariance ( void  ) const
inline
double reco::TrackBase::covariance ( int  i,
int  j 
) const
inline

(i,j)-th element of covarianve matrix ( i, j = 0, ... 4 )

Definition at line 187 of file TrackBase.h.

References covariance_, and covIndex().

187 { return covariance_[ covIndex( i, j ) ]; }
int i
Definition: DBlmapReader.cc:9
static index covIndex(index i, index j)
covariance matrix index in array
Definition: TrackBase.h:327
int j
Definition: DBlmapReader.cc:9
float covariance_[covarianceSize]
perigee 5x5 covariance matrix
Definition: TrackBase.h:304
TrackBase::index reco::TrackBase::covIndex ( index  i,
index  j 
)
inlinestatic

covariance matrix index in array

Definition at line 327 of file TrackBase.h.

References a, b, and j.

Referenced by covariance(), reco::GsfTrack::covarianceMode(), error(), and reco::GsfTrack::errorMode().

327  {
328  int a = ( i <= j ? i : j ), b = ( i <= j ? j : i );
329  return b * ( b + 1 ) / 2 + a;
330  }
int i
Definition: DBlmapReader.cc:9
int j
Definition: DBlmapReader.cc:9
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
double reco::TrackBase::d0 ( ) const
inline
double reco::TrackBase::d0Error ( ) const
inline
double reco::TrackBase::dsz ( ) const
inline

dsz parameter (THIS IS NOT the SZ impact parameter to (0,0,0) if refPoint is far from (0,0,0): see parametrization definition above for details)

Definition at line 125 of file TrackBase.h.

References p(), pt(), px(), py(), pz(), vx(), vy(), and vz().

Referenced by TestOutliers::analyze(), CosmicSplitterValidation::is_gold_muon(), MatcherByPullsAlgorithm::match(), muonisolation::CaloExtractor::MuonAtCaloPosition(), RecoTrackSelector::operator()(), parameters(), ShallowTracksProducer::produce(), RecoTracktoTP::r_dsz(), TPtoRecoTrack::rA_dsz(), and TPtoRecoTrack::rB_dsz().

125 { return vz()*pt()/p() - (vx()*px()+vy()*py())/pt() * pz()/p(); }
double p() const
momentum vector magnitude
Definition: TrackBase.h:129
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:133
double pt() const
track transverse momentum
Definition: TrackBase.h:131
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:137
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:147
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:145
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:135
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:143
double reco::TrackBase::dsz ( const Point myBeamSpot) const
inline

dsz parameter with respect to a user-given beamSpot (WARNING: this quantity can only be interpreted as the distance in the S-Z plane to the beamSpot, if the beam spot is reasonably close to the refPoint, since linear approximations are involved). This is a good approximation for Tracker tracks.

Definition at line 169 of file TrackBase.h.

References p(), pt(), px(), py(), pz(), vx(), vy(), and vz().

169  {
170  return (vz()-myBeamSpot.z())*pt()/p() - ((vx()-myBeamSpot.x())*px()+(vy()-myBeamSpot.y())*py())/pt() *pz()/p();
171  }
double p() const
momentum vector magnitude
Definition: TrackBase.h:129
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:133
double pt() const
track transverse momentum
Definition: TrackBase.h:131
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:137
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:147
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:145
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:135
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:143
double reco::TrackBase::dszError ( ) const
inline

error on dsz

Definition at line 213 of file TrackBase.h.

References error(), and i_dsz.

Referenced by ShallowTracksProducer::produce().

213 { return error( i_dsz ); }
double error(int i) const
error on specified element
Definition: TrackBase.h:189
double reco::TrackBase::dxy ( ) const
inline

dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close to (0,0,0): see parametrization definition above for details). See also function dxy(myBeamSpot) below.

Definition at line 121 of file TrackBase.h.

References pt(), px(), py(), vx(), and vy().

Referenced by EwkMuLumiMonitorDQM::analyze(), TrackAnalyzer::analyze(), CheckHitPattern::analyze(), TkConvValidator::analyze(), TestOutliers::analyze(), GenPurposeSkimmerData::analyze(), HLTMuonMatchAndPlot::analyze(), IsolatedTracksNxN::analyze(), MuonTrackValidator::analyze(), d0(), KalmanAlignmentTrackRefitter::debugTrackData(), muonisolation::PixelTrackExtractor::directionAtPresetRadius(), PFDisplacedVertexHelper::TracksSelector::dxy(), dxy(), AlignmentMonitorMuonSystemMap1D::event(), AlignmentMonitorMuonVsCurvature::event(), AlignmentMonitorSegmentDifferences::event(), MTVHistoProducerAlgoForTracker::fill_generic_recoTrack_histos(), MTVHistoProducerAlgoForTracker::fill_ResoAndPull_recoTrack_histos(), PrimaryVertexAnalyzer4PU::fillTrackHistos(), PrimaryVertexAnalyzer4PU::getSimEvents(), PFDisplacedVertexCandidateFinder::goodPtResolution(), spr::goodTrack(), HLTMuonTrackMassFilter::hltFilter(), PFDisplacedVertexHelper::isTrackSelected(), MatcherByPullsAlgorithm::match(), muonisolation::CaloExtractor::MuonAtCaloPosition(), HIPixelTrackFilter::operator()(), HIProtoTrackFilter::operator()(), reco::tau::RecoTauEnergyRecoveryPlugin::operator()(), reco::tau::RecoTauPiZeroStripPlugin2::operator()(), RecoTrackSelector::operator()(), parameters(), PrintRecoObjects::print(), ShallowTracksProducer::produce(), TrackIPProducer::produce(), reco::CentralityProducer::produce(), TriggerMatcherToHLTDebug::produce(), RecoTracktoTP::r_d02(), RecoTracktoTP::r_dxy(), TPtoRecoTrack::rA_d02(), TPtoRecoTrack::rA_dxy(), TPtoRecoTrack::rB_d02(), TPtoRecoTrack::rB_dxy(), TrackClassifier::reconstructionInformation(), MuonAlignmentFromReference::run(), pf2pat::IPCutPFCandidateSelectorDefinition::select(), reco::modules::CosmicTrackSelector::select(), reco::modules::MultiTrackSelector::select(), HLTMuonMatchAndPlot::selectedMuons(), TrackWithVertexSelector::testVertices(), reco::modules::TrackMultiSelector::testVtx(), ConversionProducer::trackD0Cut(), and QcdUeDQM::trackSelection().

121 { return ( - vx() * py() + vy() * px() ) / pt(); }
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:133
double pt() const
track transverse momentum
Definition: TrackBase.h:131
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:145
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:135
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:143
double reco::TrackBase::dxy ( const Point myBeamSpot) const
inline

dxy parameter with respect to a user-given beamSpot (WARNING: this quantity can only be interpreted as a minimum transverse distance if beamSpot, if the beam spot is reasonably close to the refPoint, since linear approximations are involved). This is a good approximation for Tracker tracks.

Definition at line 159 of file TrackBase.h.

References pt(), px(), py(), vx(), and vy().

159  {
160  return ( - (vx()-myBeamSpot.x()) * py() + (vy()-myBeamSpot.y()) * px() ) / pt();
161  }
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:133
double pt() const
track transverse momentum
Definition: TrackBase.h:131
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:145
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:135
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:143
double reco::TrackBase::dxy ( const BeamSpot theBeamSpot) const
inline

dxy parameter with respect to the beamSpot taking into account the beamspot slopes WARNING: this quantity can only be interpreted as a minimum transverse distance if beamSpot, if the beam spot is reasonably close to the refPoint, since linear approximations are involved). This is a good approximation for Tracker tracks.

Definition at line 164 of file TrackBase.h.

References dxy(), reco::BeamSpot::position(), and vz().

164  {
165  return dxy(theBeamSpot.position(vz()));
166  }
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:147
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:121
double reco::TrackBase::dxyError ( ) const
inline
double reco::TrackBase::dz ( ) const
inline

dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to (0,0,0). See also function dz(myBeamSpot) below.

Definition at line 127 of file TrackBase.h.

References pt(), px(), py(), pz(), vx(), vy(), and vz().

Referenced by TrackAnalyzer::analyze(), ValidationMisalignedTracker::analyze(), TrackSplittingMonitor::analyze(), CheckHitPattern::analyze(), TestOutliers::analyze(), HLTMuonMatchAndPlot::analyze(), CosmicSplitterValidation::analyze(), MuonTrackValidator::analyze(), IsolatedTracksNxN::analyze(), AlignmentTrackSelector::basicCuts(), KalmanAlignmentTrackRefitter::debugTrackData(), AlignmentMonitorGeneric::event(), MTVHistoProducerAlgoForTracker::fill_generic_recoTrack_histos(), MTVHistoProducerAlgoForTracker::fill_ResoAndPull_recoTrack_histos(), MillePedeMonitor::fillTrack(), TrackerValidationVariables::fillTrackQuantities(), PFTauVertexSelector::filter(), ConversionFinder::getConversionInfo(), spr::goodTrack(), HLTMuonTrackMassFilter::hltFilter(), muonisolation::CaloExtractor::MuonAtCaloPosition(), HIPixelTrackFilter::operator()(), reco::tau::RecoTauEnergyRecoveryPlugin::operator()(), reco::tau::RecoTauPiZeroStripPlugin2::operator()(), JetTracksAssociationDRVertexAssigned::produce(), TrackIPProducer::produce(), reco::CentralityProducer::produce(), TriggerMatcherToHLTDebug::produce(), fireworks::pushPixelHits(), RecoTracktoTP::r_dz(), RecoTracktoTP::r_dz2(), TPtoRecoTrack::rA_dz(), TPtoRecoTrack::rA_dz2(), TPtoRecoTrack::rB_dz(), TPtoRecoTrack::rB_dz2(), TrackClassifier::reconstructionInformation(), HIPAlignmentAlgorithm::run(), pf2pat::IPCutPFCandidateSelectorDefinition::select(), reco::modules::CosmicTrackSelector::select(), reco::modules::MultiTrackSelector::select(), HLTMuonMatchAndPlot::selectedMuons(), PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::selectPriVtxCompatibleWithTrack(), TrackWithVertexSelector::testTrack(), TrackWithVertexSelector::testVertices(), reco::modules::TrackMultiSelector::testVtx(), and QcdUeDQM::trackSelection().

127 { return vz() - (vx()*px()+vy()*py())/pt() * (pz()/pt()); }
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:133
double pt() const
track transverse momentum
Definition: TrackBase.h:131
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:137
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:147
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:145
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:135
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:143
double reco::TrackBase::dz ( const Point myBeamSpot) const
inline

dz parameter with respect to a user-given beamSpot (WARNING: this quantity can only be interpreted as the track z0, if the beamSpot is reasonably close to the refPoint, since linear approximations are involved). This is a good approximation for Tracker tracks.

Definition at line 173 of file TrackBase.h.

References pt(), px(), py(), pz(), vx(), vy(), and vz().

173  {
174  return (vz()-myBeamSpot.z()) - ((vx()-myBeamSpot.x())*px()+(vy()-myBeamSpot.y())*py())/pt() * pz()/pt();
175  }
double px() const
x coordinate of momentum vector
Definition: TrackBase.h:133
double pt() const
track transverse momentum
Definition: TrackBase.h:131
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:137
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:147
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:145
double py() const
y coordinate of momentum vector
Definition: TrackBase.h:135
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:143
double reco::TrackBase::dzError ( ) const
inline
double reco::TrackBase::error ( int  i) const
inline
double reco::TrackBase::eta ( ) const
inline

pseudorapidity of momentum vector

Definition at line 141 of file TrackBase.h.

References momentum_.

Referenced by DeDxDiscriminatorLearner::algoAnalyze(), SiStripGainFromData::algoAnalyze(), EwkMuLumiMonitorDQM::analyze(), TrackAnalyzer::analyze(), HLTMonBTagMuSource::analyze(), SegmentTrackAnalyzer::analyze(), ValidationMisalignedTracker::analyze(), CheckHitPattern::analyze(), GlobalMuonMatchAnalyzer::analyze(), TestTrackHits::analyze(), TestOutliers::analyze(), GetTrackTrajInfo::analyze(), GenPurposeSkimmerData::analyze(), CosmicSplitterValidation::analyze(), SiPixelLorentzAngle::analyze(), IsolatedTracksCone::analyze(), MuonTrackValidator::analyze(), IsolatedTracksNxN::analyze(), PrimaryVertexAnalyzer4PU::analyzeVertexCollection(), CalibrationTrackSelector::basicCuts(), AlignmentTrackSelector::basicCuts(), FWTrackProxyBuilder::build(), pat::CaloIsolationEnergy::calculate(), JetPlusTrackProducerAA::calculateBGtracksJet(), spr::chargeIsolationHcal(), compEcalEnergySum(), compHcalEnergySum(), KalmanAlignmentTrackRefitter::debugTrackData(), muonisolation::ExtractorFromDeposits::deposit(), muonisolation::CaloExtractor::deposit(), muonisolation::PixelTrackExtractor::deposit(), muonisolation::TrackExtractor::deposit(), muonisolation::JetExtractor::deposit(), muonisolation::CaloExtractorByAssociator::deposit(), PFCandWithSuperClusterExtractor::depositFromObject(), muonisolation::CaloExtractorByAssociator::deposits(), muonisolation::PixelTrackExtractor::directionAtPresetRadius(), TrackAnalyzer::doTrackerSpecificFillHists(), PrimaryVertexAnalyzer4PU::dumpHitInfo(), MuonCaloCompatibility::evaluate(), AlignmentMonitorGeneric::event(), MTVHistoProducerAlgoForTracker::fill_generic_recoTrack_histos(), MTVHistoProducerAlgoForTracker::fill_ResoAndPull_recoTrack_histos(), MTVHistoProducerAlgoForTracker::fill_simAssociated_recoTrack_histos(), TrackAnalyzer::fillHistosForState(), JPTJetAnalyzer::fillTrackHistograms(), PrimaryVertexAnalyzer4PU::fillTrackHistos(), TrackerValidationVariables::fillTrackQuantities(), reco::HcalNoiseInfoProducer::filltracks(), MuonAlignmentPreFilter::filter(), CSCEfficiency::filter(), PFElecTkProducer::FindPfRef(), MuonShowerInformationFiller::getCompatibleDets(), Tau3MuReco::getDeltaR(), MuonIdProducer::isGoodTrack(), PFDisplacedVertexCandidateFinder::link(), MatcherByPullsAlgorithm::match(), HSCPTrackSelector::matchingMuon(), TrackClassFilter::operator()(), reco::tau::RecoTauPiZeroStripPlugin2::operator()(), RecoTrackSelector::operator()(), PrimaryVertexAnalyzer4PU::printEventSummary(), IsolatedTracksCone::printTrack(), IsolatedTracksNxN::printTrack(), examples::TrackAnalysisAlgorithm::process(), ShallowTracksProducer::produce(), JetTracksAssociationDRVertex::produce(), QuarkoniaTrackSelector::produce(), JetTracksAssociationDRVertexAssigned::produce(), JetVetoedTracksAssociationDRVertex::produce(), L3MuonIsolationProducer::produce(), L3MuonCombinedRelativeIsolationProducer::produce(), reco::modules::AnalyticalTrackSelector::produce(), DeDxDiscriminatorProducer::produce(), reco::CentralityProducer::produce(), TriggerMatcherToHLTDebug::produce(), fireworks::pushPixelHits(), SoftLepton::refineJetAxis(), CutsIsolatorWithCorrection::result(), SimpleCutsIsolator::result(), MuonMillepedeAlgorithm::run(), HIPAlignmentAlgorithm::run(), MuonDTLocalMillepedeAlgorithm::run(), reco::modules::CosmicTrackSelector::select(), reco::modules::MultiTrackSelector::select(), FWTrackHitsDetailView::setTextInfo(), FWConvTrackHitsDetailView::setTextInfo(), FWPFTrackUtils::setupLegoTrack(), FWPFTrackUtils::setupTrack(), reco::SoftLeptonTagInfo::taggingVariables(), TrackWithVertexSelector::testTrack(), EwkMuLumiMonitorDQM::tkIso(), QcdUeDQM::trackSelection(), muonisolation::TrackExtractor::vetos(), muonisolation::PixelTrackExtractor::vetos(), and egammaisolation::EgammaTrackExtractor::vetos().

141 { return momentum_.Eta(); }
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:310
double reco::TrackBase::etaError ( ) const
inline

error on eta

Definition at line 205 of file TrackBase.h.

References error(), i_lambda, p(), and pt().

Referenced by CosmicSplitterValidation::analyze(), TrackAnalyzer::fillHistosForState(), and ShallowTracksProducer::produce().

205 { return error( i_lambda ) * p()/pt(); }
double p() const
momentum vector magnitude
Definition: TrackBase.h:129
double pt() const
track transverse momentum
Definition: TrackBase.h:131
double error(int i) const
error on specified element
Definition: TrackBase.h:189
TrackBase::CovarianceMatrix & TrackBase::fill ( CovarianceMatrix v) const

fill SMatrix

Definition at line 42 of file TrackBase.cc.

References covariance_, and reco::fillCovariance().

Referenced by covariance().

42  {
43  return fillCovariance( v, covariance_ );
44 }
PerigeeCovarianceMatrix & fillCovariance(PerigeeCovarianceMatrix &v, const float *data)
float covariance_[covarianceSize]
perigee 5x5 covariance matrix
Definition: TrackBase.h:304
const HitPattern& reco::TrackBase::hitPattern ( ) const
inline

Access the hit pattern, indicating in which Tracker layers the track has hits.

Definition at line 223 of file TrackBase.h.

References hitPattern_.

Referenced by TrackAnalyzer::analyze(), CheckHitPattern::analyze(), PFCheckHitPattern::analyze(), GetTrackTrajInfo::analyze(), IsolatedTracksCone::analyze(), IsolatedTracksNxN::analyze(), MuonTrackValidator::analyze(), spr::coneChargeIsolation(), TrackAnalyzer::doTrackerSpecificFillHists(), MTVHistoProducerAlgoForTracker::fill_simAssociated_recoTrack_histos(), TrackAnalyzer::fillHistosForState(), PrimaryVertexAnalyzer4PU::fillTrackHistos(), spr::goodTrack(), reco::TransientTrack::hitPattern(), HLTTrackWithHits::hltFilter(), PFMuonAlgo::isTrackerTightMuon(), reco::TransientTrack::numberOfLostHits(), reco::TransientTrack::numberOfValidHits(), reco::TrackSelector::operator()(), GhostTrackComputer::operator()(), TrackClassFilter::operator()(), reco::tau::RecoTauEnergyRecoveryPlugin::operator()(), reco::tau::RecoTauPiZeroStripPlugin2::operator()(), RecoTrackSelector::operator()(), FWTrackResidualDetailView::prepareData(), CheckHitPattern::print(), PFCheckHitPattern::print(), PFMuonAlgo::printMuonProperties(), IsolatedTracksCone::printTrack(), IsolatedTracksNxN::printTrack(), TrackIPProducer::produce(), TrackClusterRemover::produce(), PFAlgo::reconstructTrack(), reco::Track::residualX(), reco::Track::residualY(), HIPAlignmentAlgorithm::run(), reco::modules::CosmicTrackSelector::select(), reco::modules::MultiTrackSelector::select(), reco::modules::TrackMultiSelector::select(), FWConvTrackHitsDetailView::setTextInfo(), TrackWithVertexSelector::testTrack(), and QcdUeDQM::trackSelection().

223 { return hitPattern_; }
HitPattern hitPattern_
hit pattern
Definition: TrackBase.h:297
bool reco::TrackBase::isLooper ( ) const
inline

Definition at line 292 of file TrackBase.h.

References nLoops_.

292 { return (nLoops_>0);}
signed char nLoops_
number of loops made during the building of the trajectory of a looper particle
Definition: TrackBase.h:323
double reco::TrackBase::lambda ( ) const
inline

Lambda angle.

Definition at line 119 of file TrackBase.h.

References M_PI, and momentum_.

Referenced by TestOutliers::analyze(), PrimaryVertexAnalyzer4PU::fillTrackHistos(), MTVHistoProducerAlgoForTracker::getRecoMomentum(), MuonTrackValidator::getRecoMomentum(), and parameters().

119 { return M_PI/2 - momentum_.theta(); }
#define M_PI
Definition: BFit3D.cc:3
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:310
double reco::TrackBase::lambdaError ( ) const
inline

error on lambda

Definition at line 203 of file TrackBase.h.

References error(), and i_lambda.

Referenced by MTVHistoProducerAlgoForTracker::getRecoMomentum(), MuonTrackValidator::getRecoMomentum(), and PixelTrackFilterByKinematics::operator()().

203 { return error( i_lambda ); }
double error(int i) const
error on specified element
Definition: TrackBase.h:189
const Vector& reco::TrackBase::momentum ( ) const
inline

track momentum vector

Definition at line 150 of file TrackBase.h.

References momentum_.

Referenced by ValidationMisalignedTracker::analyze(), PhotonValidator::analyze(), L25TauAnalyzer::analyze(), IsolatedTracksCone::analyze(), IsolatedTracksHcalScale::analyze(), MuonTrackValidator::analyze(), IsolatedTracksNxN::analyze(), converter::StandAloneMuonTrackToCandidate::convert(), converter::TrackToCandidate::convert(), SiStripFineDelayHit::detId(), EgammaHLTTrackIsolation::electronIsolation(), MTVHistoProducerAlgoForTracker::fill_generic_recoTrack_histos(), MillePedeMonitor::fillTrack(), TrackDetectorAssociator::getFreeTrajectoryState(), HTrackAssociator::getFreeTrajectoryState(), VZeroFinder::getGlobalTrajectoryParameters(), PFDisplacedVertexCandidateFinder::getGlobalTrajectoryParameters(), ValidHitPairFilter::getTrajectory(), SeedFromProtoTrack::init(), trajectoryStateTransform::initialFreeState(), HLTmmkkFilter::initialFreeState(), HLTmmkFilter::initialFreeState(), reco::modules::CosmicTrackSplitter::makeCandidate(), reco::modules::TrackerTrackHitFilter::makeCandidate(), MuonErrorMatrixAdjuster::makeTrack(), reco::PFDisplacedVertex::momentum(), reco::TrackSelector::operator()(), CombinedSVComputer::operator()(), GhostTrackComputer::operator()(), reco::print(), PrintRecoObjects::print(), IsolatedTracksCone::printTrack(), IsolatedTracksNxN::printTrack(), SeedGeneratorFromProtoTracksEDProducer::produce(), SoftLepton::refineJetAxis(), and ObjectValidator::validTrack().

150 { return momentum_; }
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:310
double reco::TrackBase::ndof ( ) const
inline
signed char reco::TrackBase::nLoops ( ) const
inline

Definition at line 293 of file TrackBase.h.

References nLoops_.

293 {return nLoops_;}
signed char nLoops_
number of loops made during the building of the trajectory of a looper particle
Definition: TrackBase.h:323
double reco::TrackBase::normalizedChi2 ( ) const
inline

chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)

Definition at line 111 of file TrackBase.h.

References chi2_, and ndof_.

Referenced by TrackAnalyzer::analyze(), HLTMonBTagMuSource::analyze(), ValidationMisalignedTracker::analyze(), TrackSplittingMonitor::analyze(), TestOutliers::analyze(), CosmicSplitterValidation::analyze(), MuonTrackValidator::analyze(), IsolatedTracksNxN::analyze(), helper::SimpleJetTrackAssociator::associate(), helper::SimpleJetTrackAssociator::associateTransient(), CalibrationTrackSelector::basicCuts(), AlignmentTrackSelector::basicCuts(), KalmanAlignmentTrackRefitter::debugTrackData(), AlignmentMonitorGeneric::event(), AlignmentMonitorTracksFromTrajectories::event(), MTVHistoProducerAlgoForTracker::fill_ResoAndPull_recoTrack_histos(), MTVHistoProducerAlgoForTracker::fill_simAssociated_recoTrack_histos(), TrackAnalyzer::fillHistosForState(), MillePedeMonitor::fillTrack(), PrimaryVertexAnalyzer4PU::fillTrackHistos(), TrackerValidationVariables::fillTrackQuantities(), CSCEfficiency::filter(), PFDisplacedVertexCandidateFinder::goodPtResolution(), spr::goodTrack(), HLTMuonTrackMassFilter::hltFilter(), PFDisplacedVertexHelper::isTrackSelected(), MuonResidualsFromTrack::normalizedChi2(), reco::TrackSelector::operator()(), GhostTrackComputer::operator()(), TrackClassFilter::operator()(), reco::tau::RecoTauEnergyRecoveryPlugin::operator()(), reco::tau::RecoTauPiZeroStripPlugin2::operator()(), RecoTrackSelector::operator()(), IsolatedTracksCone::printTrack(), IsolatedTracksNxN::printTrack(), AlignmentMonitorMuonVsCurvature::processMuonResidualsFromTrack(), TrackIPProducer::produce(), MuonMillepedeAlgorithm::run(), HIPAlignmentAlgorithm::run(), reco::modules::CosmicTrackSelector::select(), reco::modules::MultiTrackSelector::select(), reco::modules::TrackMultiSelector::select(), reco::SoftLeptonTagInfo::taggingVariables(), TrackWithVertexSelector::testTrack(), ConversionProducer::trackQualityFilter(), and QcdUeDQM::trackSelection().

111 { return ndof_ != 0 ? chi2_ / ndof_ : chi2_ * 1e6; }
float chi2_
chi-squared
Definition: TrackBase.h:306
float ndof_
number of degrees of freedom
Definition: TrackBase.h:313
unsigned short reco::TrackBase::numberOfLostHits ( ) const
inline
unsigned short reco::TrackBase::numberOfValidHits ( ) const
inline

number of valid hits found

Definition at line 232 of file TrackBase.h.

References hitPattern_, and reco::HitPattern::numberOfValidHits().

Referenced by TrackAnalyzer::analyze(), HLTMonBTagMuSource::analyze(), ValidationMisalignedTracker::analyze(), TestTrackHits::analyze(), TestOutliers::analyze(), MuonTrackValidator::analyze(), helper::SimpleJetTrackAssociator::associate(), helper::SimpleJetTrackAssociator::associateTransient(), CalibrationTrackSelector::basicCuts(), AlignmentTrackSelector::basicCuts(), MuonAlignmentAnalyzer::doMatching(), MTVHistoProducerAlgoForTracker::fill_generic_recoTrack_histos(), MTVHistoProducerAlgoForTracker::fill_recoAssociated_simTrack_histos(), MTVHistoProducerAlgoForTracker::fill_ResoAndPull_recoTrack_histos(), MTVHistoProducerAlgoForTracker::fill_simAssociated_recoTrack_histos(), MillePedeMonitor::fillTrack(), TrackerValidationVariables::fillTrackQuantities(), reco::Track::found(), HLTMuonTrackMassFilter::hltFilter(), PFDisplacedVertexHelper::isTrackSelected(), TrackClassFilter::operator()(), IsolatedTracksCone::printTrack(), IsolatedTracksNxN::printTrack(), ShallowTracksProducer::produce(), cms::TrackListMerger::produce(), TriggerMatcherToHLTDebug::produce(), MuonMillepedeAlgorithm::run(), HIPAlignmentAlgorithm::run(), reco::modules::MultiTrackSelector::select(), TrackWithVertexSelector::testTrack(), and ObjectValidator::validTrack().

232 { return hitPattern_.numberOfValidHits(); }
int numberOfValidHits() const
Definition: HitPattern.h:554
HitPattern hitPattern_
hit pattern
Definition: TrackBase.h:297
double reco::TrackBase::p ( ) const
inline

momentum vector magnitude

Definition at line 129 of file TrackBase.h.

References momentum_.

Referenced by PFTrackTransformer::addPoints(), DeDxDiscriminatorLearner::algoAnalyze(), SiStripGainFromData::algoAnalyze(), EwkMuLumiMonitorDQM::analyze(), CosmicSplitterValidation::analyze(), IsolatedTracksCone::analyze(), IsolatedTracksHcalScale::analyze(), IsolatedTracksNxN::analyze(), PrimaryVertexAnalyzer4PU::analyzeVertexCollection(), AlignmentTrackSelector::basicCuts(), HBHEHitMap::calcTracksNeighborTowers_(), HBHEHitMap::calcTracksSameTowers_(), pat::TrackerIsolationPt::calculate(), spr::chargeIsolation(), spr::chargeIsolationCone(), spr::chargeIsolationEcal(), spr::chargeIsolationHcal(), spr::coneChargeIsolation(), PFRecoTauDiscriminationAgainstMuon2::discriminate(), dsz(), dzError(), etaError(), MuonCaloCompatibility::evaluate(), AlignmentMonitorMuonVsCurvature::event(), AlignmentMonitorMuonSystemMap1D::event(), AlignmentMonitorSegmentDifferences::event(), TrackAnalyzer::fillHistosForState(), TrackerValidationVariables::fillTrackQuantities(), reco::HcalNoiseInfoProducer::filltracks(), MuonAlignmentPreFilter::filter(), CSCEfficiency::filter(), MuonShowerInformationFiller::getCompatibleDets(), ConversionFinder::getConversionInfo(), MuonIdProducer::isGoodTrack(), PFDisplacedVertexHelper::lambdaCP(), MuonIdProducer::makeMuon(), MuonIdTruthInfo::matchChi2(), muonid::matchTracks(), muonisolation::CaloExtractor::MuonAtCaloPosition(), TrackClassFilter::operator()(), reco::tau::RecoTauEnergyRecoveryPlugin::operator()(), IsolatedTracksCone::printTrack(), IsolatedTracksNxN::printTrack(), ShallowTracksProducer::produce(), JetTracksAssociationDRVertex::produce(), QuarkoniaTrackSelector::produce(), JetTracksAssociationDRVertexAssigned::produce(), JetVetoedTracksAssociationDRVertex::produce(), SiStripElectronAssociator::produce(), DeDxDiscriminatorProducer::produce(), ptError(), qoverp(), PFAlgo::reconstructTrack(), MuonTrackingRegionBuilder::region(), HIPAlignmentAlgorithm::run(), MuonDTLocalMillepedeAlgorithm::run(), CSCOverlapsAlignmentAlgorithm::run(), MuonAlignmentFromReference::run(), CandCommonVertexFitterBase::set(), PFCandCommonVertexFitterBase::set(), reco::SoftLeptonTagInfo::taggingVariables(), and reco::JetTracksAssociation::tracksP4().

129 { return momentum_.R(); }
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:310
double reco::TrackBase::parameter ( int  i) const
inline

i-th parameter ( i = 0, ... 4 )

Definition at line 185 of file TrackBase.h.

References i, and parameters().

Referenced by TriggerMatcherToHLTDebug::produce().

185 { return parameters()[i]; }
int i
Definition: DBlmapReader.cc:9
ParameterVector parameters() const
Track parameters with one-to-one correspondence to the covariance matrix.
Definition: TrackBase.h:178
ParameterVector reco::TrackBase::parameters ( void  ) const
inline

Track parameters with one-to-one correspondence to the covariance matrix.

Definition at line 178 of file TrackBase.h.

References dsz(), dxy(), lambda(), phi(), and qoverp().

Referenced by MuonTrackValidator::analyze(), SiStripFineDelayHit::detId(), Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::inputTags(), parameter(), Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::properties(), and Vispa.Plugins.ConfigEditor.ConfigDataAccessor.ConfigDataAccessor::recursePSetProperties().

178  {
179  return ParameterVector(qoverp(),lambda(),phi(),dxy(),dsz());
180  }
double qoverp() const
q/p
Definition: TrackBase.h:115
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:139
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:70
double dsz() const
dsz parameter (THIS IS NOT the SZ impact parameter to (0,0,0) if refPoint is far from (0...
Definition: TrackBase.h:125
double lambda() const
Lambda angle.
Definition: TrackBase.h:119
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:121
double reco::TrackBase::phi ( ) const
inline

azimuthal angle of momentum vector

Definition at line 139 of file TrackBase.h.

References momentum_.

Referenced by TrackAnalyzer::analyze(), HLTMonBTagMuSource::analyze(), ValidationMisalignedTracker::analyze(), SegmentTrackAnalyzer::analyze(), TrackSplittingMonitor::analyze(), TestOutliers::analyze(), GenPurposeSkimmerData::analyze(), CosmicSplitterValidation::analyze(), SiPixelLorentzAngle::analyze(), PrimaryVertexAnalyzer4PU::analyzeVertexCollection(), CalibrationTrackSelector::basicCuts(), AlignmentTrackSelector::basicCuts(), pat::CaloIsolationEnergy::calculate(), JetPlusTrackProducerAA::calculateBGtracksJet(), spr::chargeIsolationHcal(), AlignmentTwoBodyDecayTrackSelector::checkAcoplanarity(), AlignmentTwoBodyDecayTrackSelector::checkMETAcoplanarity(), compEcalEnergySum(), compHcalEnergySum(), KalmanAlignmentTrackRefitter::debugTrackData(), muonisolation::ExtractorFromDeposits::deposit(), muonisolation::CaloExtractor::deposit(), muonisolation::PixelTrackExtractor::deposit(), muonisolation::TrackExtractor::deposit(), muonisolation::JetExtractor::deposit(), muonisolation::CaloExtractorByAssociator::deposit(), PFCandWithSuperClusterExtractor::depositFromObject(), muonisolation::CaloExtractorByAssociator::deposits(), muonisolation::PixelTrackExtractor::directionAtPresetRadius(), TrackAnalyzer::doTrackerSpecificFillHists(), AlignmentMonitorGeneric::event(), TrackAnalyzer::fillHistosForState(), JPTJetAnalyzer::fillTrackHistograms(), PrimaryVertexAnalyzer4PU::fillTrackHistos(), TrackerValidationVariables::fillTrackQuantities(), CSCEfficiency::filter(), PFElecTkProducer::FindPfRef(), MuonShowerInformationFiller::getCompatibleDets(), Tau3MuReco::getDeltaR(), MTVHistoProducerAlgoForTracker::getRecoMomentum(), MuonTrackValidator::getRecoMomentum(), CosmicSplitterValidation::is_gold_muon(), PFDisplacedVertexCandidateFinder::link(), MatcherByPullsAlgorithm::match(), MuonIdTruthInfo::matchChi2(), HSCPTrackSelector::matchingMuon(), muonisolation::CaloExtractor::MuonAtCaloPosition(), reco::tau::RecoTauPiZeroStripPlugin2::operator()(), parameters(), PrimaryVertexAnalyzer4PU::printEventSummary(), IsolatedTracksCone::printTrack(), IsolatedTracksNxN::printTrack(), ShallowTracksProducer::produce(), JetTracksAssociationDRVertex::produce(), JetTracksAssociationDRVertexAssigned::produce(), JetVetoedTracksAssociationDRVertex::produce(), SoftLepton::refineJetAxis(), MuonMillepedeAlgorithm::run(), HIPAlignmentAlgorithm::run(), MuonDTLocalMillepedeAlgorithm::run(), FWTrackHitsDetailView::setTextInfo(), FWConvTrackHitsDetailView::setTextInfo(), reco::SoftLeptonTagInfo::taggingVariables(), EwkMuLumiMonitorDQM::tkIso(), QcdUeDQM::trackSelection(), muonisolation::TrackExtractor::vetos(), muonisolation::PixelTrackExtractor::vetos(), and egammaisolation::EgammaTrackExtractor::vetos().

139 { return momentum_.Phi(); }
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:310
double reco::TrackBase::phiError ( ) const
inline
double reco::TrackBase::pt ( ) const
inline

track transverse momentum

Definition at line 131 of file TrackBase.h.

References momentum_, and mathSSE::sqrt().

Referenced by PFTrackTransformer::addPoints(), PFTrackTransformer::addPointsAndBrems(), SiStripGainFromData::algoAnalyze(), PFCandConnector::analyseNuclearWSec(), MultiTrackValidator::analyze(), EwkMuLumiMonitorDQM::analyze(), TrackAnalyzer::analyze(), HLTMonBTagMuSource::analyze(), ValidationMisalignedTracker::analyze(), SegmentTrackAnalyzer::analyze(), TrackSplittingMonitor::analyze(), GlobalMuonMatchAnalyzer::analyze(), ElectronConversionRejectionValidator::analyze(), TkConvValidator::analyze(), TestTrackHits::analyze(), L25TauAnalyzer::analyze(), TestOutliers::analyze(), GenPurposeSkimmerData::analyze(), CosmicSplitterValidation::analyze(), SiPixelLorentzAngle::analyze(), IsolatedTracksCone::analyze(), IsolatedTracksHcalScale::analyze(), MuonTrackValidator::analyze(), IsolatedTracksNxN::analyze(), SoftLeptonTagPlotter::analyzeTag(), CalibrationTrackSelector::basicCuts(), AlignmentTrackSelector::basicCuts(), FWTrackProxyBuilder::build(), TrackProducerAlgorithm< reco::Track >::buildTrack(), Tau3MuReco::check4MuonTrack(), KalmanAlignmentTrackRefitter::debugTrackData(), muonisolation::CaloExtractor::deposit(), muonisolation::PixelTrackExtractor::deposit(), muonisolation::TrackExtractor::deposit(), PFCandWithSuperClusterExtractor::depositFromObject(), muonisolation::PixelTrackExtractor::directionAtPresetRadius(), dsz(), PrimaryVertexAnalyzer4PU::dumpHitInfo(), dxy(), dz(), dzError(), fw::estimate_field(), etaError(), AlignmentMonitorMuonVsCurvature::event(), AlignmentMonitorMuonSystemMap1D::event(), AlignmentMonitorSegmentDifferences::event(), AlignmentMonitorGeneric::event(), TrackAnalyzer::fillHistosForState(), MillePedeMonitor::fillTrack(), JPTJetAnalyzer::fillTrackHistograms(), PrimaryVertexAnalyzer4PU::fillTrackHistos(), TrackerValidationVariables::fillTrackQuantities(), reco::HcalNoiseInfoProducer::filltracks(), MuonAlignmentPreFilter::filter(), CSCEfficiency::filter(), ConversionFinder::getConversionInfo(), MTVHistoProducerAlgoForTracker::getRecoMomentum(), MuonTrackValidator::getRecoMomentum(), JPTJetTesterUnCorr::getSumPt(), JPTJetTester::getSumPt(), PFDisplacedVertexCandidateFinder::goodPtResolution(), spr::goodTrack(), MuonIdProducer::isGoodTrack(), PFDisplacedVertexHelper::isTrackSelected(), PFDisplacedVertexCandidateFinder::link(), MatcherByPullsAlgorithm::match(), HSCPTrackSelector::matchingMuon(), muonid::matchTracks(), muonisolation::CaloExtractor::MuonAtCaloPosition(), PixelTrackFilterByKinematics::operator()(), HIPixelTrackFilter::operator()(), HIProtoTrackFilter::operator()(), reco::TrackSelector::operator()(), reco::tau::RecoTauEnergyRecoveryPlugin::operator()(), CalibrationTrackSelector::ComparePt::operator()(), reco::tau::RecoTauPiZeroStripPlugin2::operator()(), AlignmentTrackSelector::ComparePt::operator()(), RecoTrackSelector::operator()(), reco::Vertex::TrackEqual::operator()(), fireworks::prepareTrack(), PrimaryVertexAnalyzer4PU::printEventSummary(), IsolatedTracksCone::printTrack(), IsolatedTracksNxN::printTrack(), examples::TrackAnalysisAlgorithm::process(), AlignmentMonitorSegmentDifferences::processMuonResidualsFromTrack(), AlignmentMonitorMuonVsCurvature::processMuonResidualsFromTrack(), AlignmentMonitorMuonSystemMap1D::processMuonResidualsFromTrack(), MuonAlignmentFromReference::processMuonResidualsFromTrack(), ShallowTracksProducer::produce(), QuarkoniaTrackSelector::produce(), TrackIPProducer::produce(), L3MuonCombinedRelativeIsolationProducer::produce(), reco::modules::MultiTrackSelector::produce(), reco::modules::AnalyticalTrackSelector::produce(), reco::CentralityProducer::produce(), TriggerMatcherToHLTDebug::produce(), ptError(), SoftLepton::refineJetAxis(), GlobalMuonRefitter::refit(), CutsIsolatorWithCorrection::result(), MuonMillepedeAlgorithm::run(), HIPAlignmentAlgorithm::run(), MuonDTLocalMillepedeAlgorithm::run(), MuonAlignmentFromReference::run(), reco::modules::CosmicTrackSelector::select(), reco::modules::MultiTrackSelector::select(), reco::modules::TrackMultiSelector::select(), reco::modules::HICaloCompatibleTrackSelector::selectByPFCands(), reco::modules::HICaloCompatibleTrackSelector::selectByTowers(), PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::selectPriVtxCompatibleWithTrack(), FWTrackHitsDetailView::setTextInfo(), FWConvTrackHitsDetailView::setTextInfo(), FWPFTrackUtils::setupLegoTrack(), FWPFTrackUtils::setupTrack(), TrackClusterSplitter::splitCluster(), TrackWithVertexSelector::testTrack(), EwkMuLumiMonitorDQM::tkIso(), reco::PFDisplacedVertex::trackPosition(), QcdUeDQM::trackSelection(), and ObjectValidator::validTrack().

131 { return sqrt( momentum_.Perp2() ); }
T sqrt(T t)
Definition: SSEVec.h:46
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:310
double reco::TrackBase::ptError ( ) const
inline

error on Pt (set to 1000 TeV if charge==0 for safety)

Definition at line 194 of file TrackBase.h.

References charge(), covariance(), i_lambda, i_qoverp, p(), pt(), pz(), and mathSSE::sqrt().

Referenced by TrackSplittingMonitor::analyze(), TestOutliers::analyze(), CosmicSplitterValidation::analyze(), TrackAnalyzer::fillHistosForState(), TrackerValidationVariables::fillTrackQuantities(), CSCEfficiency::filter(), MTVHistoProducerAlgoForTracker::getRecoMomentum(), MuonTrackValidator::getRecoMomentum(), PFDisplacedVertexCandidateFinder::goodPtResolution(), reco::tau::RecoTauEnergyRecoveryPlugin::operator()(), reco::tau::RecoTauPiZeroStripPlugin2::operator()(), ShallowTracksProducer::produce(), reco::CentralityProducer::produce(), PFAlgo::reconstructTrack(), reco::modules::MultiTrackSelector::select(), TrackWithVertexSelector::testTrack(), and QcdUeDQM::trackSelection().

194  {
195  return (charge()!=0) ? sqrt(
196  pt()*pt()*p()*p()/charge()/charge() * covariance(i_qoverp,i_qoverp)
197  + 2*pt()*p()/charge()*pz() * covariance(i_qoverp,i_lambda)
198  + pz()*pz() * covariance(i_lambda,i_lambda) ) : 1.e6;
199  }
double p() const
momentum vector magnitude
Definition: TrackBase.h:129
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:182
T sqrt(T t)
Definition: SSEVec.h:46
double pt() const
track transverse momentum
Definition: TrackBase.h:131
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:137
int charge() const
track electric charge
Definition: TrackBase.h:113
double reco::TrackBase::px ( ) const
inline

x coordinate of momentum vector

Definition at line 133 of file TrackBase.h.

References momentum_.

Referenced by reco::TrackKinematics::add(), PFTrackTransformer::addPoints(), EwkMuLumiMonitorDQM::analyze(), PrimaryVertexAnalyzer4PU::analyzeVertexCollection(), FWSecVertexProxyBuilder::build(), FWVertexProxyBuilder::build(), FWTrackProxyBuilderFF::buildTrack(), pat::TrackerIsolationPt::calculate(), FWInvMassDialog::Calculate(), compEcalEnergySum(), compHcalEnergySum(), L6SLBCorrector::correction(), directionAlongMomentum(), dsz(), dummyPrediction(), dxy(), dz(), fw::estimate_field(), TrackAnalyzer::fillHistosForState(), TrackerValidationVariables::fillTrackQuantities(), ConversionFinder::getConversionInfo(), FWPFTrackUtils::getTrack(), PFDisplacedVertexHelper::lambdaCP(), MuonIdProducer::makeMuon(), muonid::matchTracks(), fireworks::prepareTrack(), QuarkoniaTrackSelector::produce(), SiStripElectronAssociator::produce(), spr::propagateECAL(), spr::propagateHCAL(), spr::propagateTracker(), spr::propagateTrackerEnd(), spr::propagateTrackToECAL(), spr::propagateTrackToHCAL(), TrackExtrapolator::propagateTrackToVolume(), PFAlgo::reconstructTrack(), PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::selectPriVtxCompatibleWithTrack(), CandCommonVertexFitterBase::set(), PFCandCommonVertexFitterBase::set(), and reco::JetTracksAssociation::tracksP4().

133 { return momentum_.x(); }
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:310
double reco::TrackBase::py ( ) const
inline

y coordinate of momentum vector

Definition at line 135 of file TrackBase.h.

References momentum_.

Referenced by reco::TrackKinematics::add(), PFTrackTransformer::addPoints(), EwkMuLumiMonitorDQM::analyze(), PrimaryVertexAnalyzer4PU::analyzeVertexCollection(), FWSecVertexProxyBuilder::build(), FWVertexProxyBuilder::build(), FWTrackProxyBuilderFF::buildTrack(), pat::TrackerIsolationPt::calculate(), FWInvMassDialog::Calculate(), compEcalEnergySum(), compHcalEnergySum(), L6SLBCorrector::correction(), directionAlongMomentum(), dsz(), dummyPrediction(), dxy(), dz(), fw::estimate_field(), TrackAnalyzer::fillHistosForState(), TrackerValidationVariables::fillTrackQuantities(), ConversionFinder::getConversionInfo(), FWPFTrackUtils::getTrack(), PFDisplacedVertexHelper::lambdaCP(), MuonIdProducer::makeMuon(), muonid::matchTracks(), fireworks::prepareTrack(), QuarkoniaTrackSelector::produce(), SiStripElectronAssociator::produce(), spr::propagateECAL(), spr::propagateHCAL(), spr::propagateTracker(), spr::propagateTrackerEnd(), spr::propagateTrackToECAL(), spr::propagateTrackToHCAL(), TrackExtrapolator::propagateTrackToVolume(), PFAlgo::reconstructTrack(), PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::selectPriVtxCompatibleWithTrack(), CandCommonVertexFitterBase::set(), PFCandCommonVertexFitterBase::set(), and reco::JetTracksAssociation::tracksP4().

135 { return momentum_.y(); }
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:310
double reco::TrackBase::pz ( ) const
inline

z coordinate of momentum vector

Definition at line 137 of file TrackBase.h.

References momentum_.

Referenced by reco::TrackKinematics::add(), PFTrackTransformer::addPoints(), EwkMuLumiMonitorDQM::analyze(), PrimaryVertexAnalyzer4PU::analyzeVertexCollection(), FWSecVertexProxyBuilder::build(), FWVertexProxyBuilder::build(), FWTrackProxyBuilderFF::buildTrack(), pat::TrackerIsolationPt::calculate(), FWInvMassDialog::Calculate(), compEcalEnergySum(), compHcalEnergySum(), L6SLBCorrector::correction(), dsz(), dummyPrediction(), dz(), TrackAnalyzer::fillHistosForState(), TrackerValidationVariables::fillTrackQuantities(), ConversionFinder::getConversionInfo(), FWPFTrackUtils::getTrack(), PFDisplacedVertexHelper::lambdaCP(), MuonIdProducer::makeMuon(), muonid::matchTracks(), fireworks::prepareTrack(), AlignmentMonitorMuonSystemMap1D::processMuonResidualsFromTrack(), AlignmentMonitorMuonVsCurvature::processMuonResidualsFromTrack(), AlignmentMonitorSegmentDifferences::processMuonResidualsFromTrack(), MuonAlignmentFromReference::processMuonResidualsFromTrack(), QuarkoniaTrackSelector::produce(), SiStripElectronAssociator::produce(), spr::propagateECAL(), spr::propagateHCAL(), spr::propagateTracker(), spr::propagateTrackerEnd(), spr::propagateTrackToECAL(), spr::propagateTrackToHCAL(), TrackExtrapolator::propagateTrackToVolume(), ptError(), PFAlgo::reconstructTrack(), PhotonConversionTrajectorySeedProducerFromSingleLegAlgo::selectPriVtxCompatibleWithTrack(), CandCommonVertexFitterBase::set(), PFCandCommonVertexFitterBase::set(), and reco::JetTracksAssociation::tracksP4().

137 { return momentum_.z(); }
Vector momentum_
momentum vector at innermost point
Definition: TrackBase.h:310
double reco::TrackBase::qoverp ( ) const
inline
double reco::TrackBase::qoverpError ( ) const
inline
bool reco::TrackBase::quality ( const TrackQuality  ) const
inline

Track quality.

Definition at line 377 of file TrackBase.h.

References confirmed, goodIterative, highPurity, lumiQueryAPI::q, quality_, and undefQuality.

Referenced by TrackAnalyzer::analyze(), IsolatedTracksCone::analyze(), IsolatedTracksNxN::analyze(), spr::chargeIsolation(), spr::chargeIsolationEcal(), spr::chargeIsolationHcal(), spr::coneChargeIsolation(), TrackAnalyzer::fillHistosForState(), CSCEfficiency::filter(), spr::goodTrack(), AlignmentTrackSelector::isOkTrkQuality(), PFDisplacedVertexHelper::isTrackSelected(), reco::TrackSelector::operator()(), TrackFilterForPVFinding::operator()(), TrackClassFilter::operator()(), RecoTrackSelector::operator()(), IsolatedTracksCone::printTrack(), IsolatedTracksNxN::printTrack(), TrackClusterRemover::produce(), reco::CentralityProducer::produce(), spr::propagateCALO(), reco::modules::HICaloCompatibleTrackSelector::selectByPFCands(), reco::modules::HICaloCompatibleTrackSelector::selectByTowers(), TrackWithVertexSelector::testTrack(), and QcdUeDQM::trackSelection().

377  {
378  switch(q){
379  case undefQuality: return (quality_==0);
380  case goodIterative:
381  {
382  return ( ((quality_ & (1<<TrackBase::confirmed))>>TrackBase::confirmed) ||
384  }
385  default:
386  {
387  return (quality_ & (1<<q))>>q;
388  }
389  }
390  return false;
391  }
uint8_t quality_
track quality
Definition: TrackBase.h:320
TrackBase::TrackQuality TrackBase::qualityByName ( const std::string &  name)
static

Definition at line 46 of file TrackBase.cc.

References spr::find(), qualityNames, qualitySize, findQualityFiles::size, and undefQuality.

Referenced by AlignmentTrackSelector::AlignmentTrackSelector(), reco::modules::AnalyticalTrackSelector::AnalyticalTrackSelector(), PrimaryVertexValidation::analyze(), IsolatedTracksCone::analyze(), IsolatedTracksNxN::analyze(), BeamFitter::BeamFitter(), CaloRecoTauTagInfoAlgorithm::CaloRecoTauTagInfoAlgorithm(), spr::chargeIsolation(), spr::chargeIsolationEcal(), spr::chargeIsolationHcal(), spr::coneChargeIsolation(), reco::modules::CosmicTrackSelector::CosmicTrackSelector(), FilterOutScraping::filter(), IsolatedTracksHcalScale::IsolatedTracksHcalScale(), PFDisplacedVertexHelper::isTrackSelected(), JetBProbabilityComputer::JetBProbabilityComputer(), JetPlusTrackProducerAA::JetPlusTrackProducerAA(), JetProbabilityComputer::JetProbabilityComputer(), LightPFTrackProducer::LightPFTrackProducer(), reco::modules::MultiTrackSelector::MultiTrackSelector(), TrackClassFilter::operator()(), PFTrackProducer::PFTrackProducer(), IsolatedTracksCone::printTrack(), IsolatedTracksNxN::printTrack(), FastTrackMerger::produce(), cms::SimpleTrackListMerger::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), PromptTrackCountingComputer::PromptTrackCountingComputer(), spr::propagateCALO(), QcdUeDQM::QcdUeDQM(), QualityCutsAnalyzer::QualityCutsAnalyzer(), RecoTrackSelector::RecoTrackSelector(), SeedClusterRemover::SeedClusterRemover(), reco::modules::HICaloCompatibleTrackSelector::selectByPFCands(), reco::modules::HICaloCompatibleTrackSelector::selectByTowers(), TrackWithVertexSelector::testTrack(), TrackClusterRemover::TrackClusterRemover(), TrackCountingComputer::TrackCountingComputer(), TrackExtrapolator::TrackExtrapolator(), TrackFilterForPVFinding::TrackFilterForPVFinding(), cms::TrackListMerger::TrackListMerger(), reco::TrackSelector::TrackSelector(), V0Fitter::V0Fitter(), and VertexAssociatorByTracks::VertexAssociatorByTracks().

46  {
49  if(index == size) return undefQuality; // better this or throw() ?
50 
51  // cast
52  return TrackQuality(index);
53 }
unsigned int index
index type
Definition: TrackBase.h:80
TrackQuality
track quality
Definition: TrackBase.h:95
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
static const std::string qualityNames[]
Definition: TrackBase.h:96
This class analyses the reconstruction quality for a given track.
Definition: TrackQuality.h:26
tuple size
Write out results.
int reco::TrackBase::qualityMask ( ) const
inline

Definition at line 286 of file TrackBase.h.

References quality_.

Referenced by PrimaryVertexAnalyzer4PU::fillTrackHistos(), and reco::modules::MultiTrackSelector::produce().

286 { return quality_; }
uint8_t quality_
track quality
Definition: TrackBase.h:320
std::string reco::TrackBase::qualityName ( TrackQuality  q)
inlinestatic

Definition at line 406 of file TrackBase.h.

References qualityNames, and qualitySize.

Referenced by TkAlCaRecoMonitor::beginJob(), IsolatedTracksCone::printTrack(), and IsolatedTracksNxN::printTrack().

406  {
407  if(int(q) < int(qualitySize) && int(q)>=0) return qualityNames[int(q)];
408  return "undefQuality";
409  }
static const std::string qualityNames[]
Definition: TrackBase.h:96
const Point& reco::TrackBase::referencePoint ( ) const
inline

Reference point on the track.

Definition at line 153 of file TrackBase.h.

References vertex_.

Referenced by TrackAnalyzer::analyze(), CosmicSplitterValidation::analyze(), SeedFromProtoTrack::init(), MuonErrorMatrixAdjuster::makeTrack(), PrintRecoObjects::print(), IsolatedTracksCone::printTrack(), and IsolatedTracksNxN::printTrack().

153 { return vertex_; }
Point vertex_
innermost (reference) point on track
Definition: TrackBase.h:308
void reco::TrackBase::setAlgorithm ( const TrackAlgorithm  a,
bool  set = true 
)
inline

position index

Track algorithm

Definition at line 273 of file TrackBase.h.

References a, algorithm_, runtimedef::set(), and undefAlgorithm.

Referenced by TrackProducerAlgorithm< reco::GsfTrack >::buildTrack(), and TrackListCombiner::produce().

double a
Definition: hdecay.h:121
uint8_t algorithm_
track algorithm
Definition: TrackBase.h:318
void set(const std::string &name, int value)
set the flag, with a run-time name
template<typename C >
void reco::TrackBase::setHitPattern ( const C &  c)
inline
template<typename I >
void reco::TrackBase::setHitPattern ( const I &  begin,
const I &  end 
)
inline

Definition at line 253 of file TrackBase.h.

References hitPattern_, and reco::HitPattern::set().

253 { hitPattern_.set( begin, end ); }
void set(const I &begin, const I &end)
Definition: HitPattern.h:139
#define end
Definition: vmac.h:38
#define begin
Definition: vmac.h:31
HitPattern hitPattern_
hit pattern
Definition: TrackBase.h:297
void reco::TrackBase::setHitPattern ( const TrackingRecHit hit,
size_t  i 
)
inline

set hit pattern for specified hit

Definition at line 260 of file TrackBase.h.

References hitPattern_, and reco::HitPattern::set().

260 { hitPattern_.set( hit, i ); }
int i
Definition: DBlmapReader.cc:9
void set(const I &begin, const I &end)
Definition: HitPattern.h:139
HitPattern hitPattern_
hit pattern
Definition: TrackBase.h:297
void reco::TrackBase::setHitPattern ( const HitPattern hitP)
inline

set hitPattern from pre-defined hitPattern

Definition at line 266 of file TrackBase.h.

References hitPattern_.

266 {hitPattern_ = hitP;}
HitPattern hitPattern_
hit pattern
Definition: TrackBase.h:297
void reco::TrackBase::setNLoops ( signed char  value)
inline

Definition at line 290 of file TrackBase.h.

References nLoops_, and relativeConstraints::value.

Referenced by TrackProducerAlgorithm< reco::Track >::buildTrack().

290 { nLoops_=value;}
signed char nLoops_
number of loops made during the building of the trajectory of a looper particle
Definition: TrackBase.h:323
void reco::TrackBase::setQuality ( const TrackQuality  ,
bool  set = true 
)
inline

Definition at line 393 of file TrackBase.h.

References lumiQueryAPI::q, quality_, and undefQuality.

Referenced by TrackBase().

393  {
394  switch(q){
395  case undefQuality: quality_=0;
396  default:
397  {
398  if (set)//regular OR if setting value to true
399  quality_= quality_ | (1<<q) ;
400  else // doing "half-XOR" if unsetting value
401  quality_= quality_ & (~(1<<q));
402  }
403  }
404  }
uint8_t quality_
track quality
Definition: TrackBase.h:320
void set(const std::string &name, int value)
set the flag, with a run-time name
void reco::TrackBase::setQualityMask ( int  qualMask)
inline

Definition at line 287 of file TrackBase.h.

References quality_.

Referenced by TrackProducerAlgorithm< reco::Track >::buildTrack().

287 {quality_ = qualMask;}
uint8_t quality_
track quality
Definition: TrackBase.h:320
template<typename C >
void reco::TrackBase::setTrackerExpectedHitsInner ( const C &  c)
inline

Definition at line 248 of file TrackBase.h.

References reco::HitPattern::set(), and trackerExpectedHitsInner_.

248 { trackerExpectedHitsInner_.set( c.begin(), c.end() ); }
HitPattern trackerExpectedHitsInner_
hit pattern used for expected crossed layers after the last trajectory&#39;s hit
Definition: TrackBase.h:299
void set(const I &begin, const I &end)
Definition: HitPattern.h:139
template<typename I >
void reco::TrackBase::setTrackerExpectedHitsInner ( const I &  begin,
const I &  end 
)
inline

Definition at line 255 of file TrackBase.h.

References reco::HitPattern::set(), and trackerExpectedHitsInner_.

HitPattern trackerExpectedHitsInner_
hit pattern used for expected crossed layers after the last trajectory&#39;s hit
Definition: TrackBase.h:299
void set(const I &begin, const I &end)
Definition: HitPattern.h:139
#define end
Definition: vmac.h:38
#define begin
Definition: vmac.h:31
void reco::TrackBase::setTrackerExpectedHitsInner ( const TrackingRecHit hit,
size_t  i 
)
inline

Definition at line 262 of file TrackBase.h.

References reco::HitPattern::set(), and trackerExpectedHitsInner_.

262 { trackerExpectedHitsInner_.set( hit, i ); }
int i
Definition: DBlmapReader.cc:9
HitPattern trackerExpectedHitsInner_
hit pattern used for expected crossed layers after the last trajectory&#39;s hit
Definition: TrackBase.h:299
void set(const I &begin, const I &end)
Definition: HitPattern.h:139
void reco::TrackBase::setTrackerExpectedHitsInner ( const HitPattern hitP)
inline

Definition at line 267 of file TrackBase.h.

References trackerExpectedHitsInner_.

HitPattern trackerExpectedHitsInner_
hit pattern used for expected crossed layers after the last trajectory&#39;s hit
Definition: TrackBase.h:299
template<typename C >
void reco::TrackBase::setTrackerExpectedHitsOuter ( const C &  c)
inline

Definition at line 250 of file TrackBase.h.

References reco::HitPattern::set(), and trackerExpectedHitsOuter_.

250 { trackerExpectedHitsOuter_.set( c.begin(), c.end() ); }
void set(const I &begin, const I &end)
Definition: HitPattern.h:139
HitPattern trackerExpectedHitsOuter_
hit pattern used for expected crossed layers before the first trajectory&#39;s hit
Definition: TrackBase.h:301
template<typename I >
void reco::TrackBase::setTrackerExpectedHitsOuter ( const I &  begin,
const I &  end 
)
inline

Definition at line 257 of file TrackBase.h.

References reco::HitPattern::set(), and trackerExpectedHitsOuter_.

void set(const I &begin, const I &end)
Definition: HitPattern.h:139
#define end
Definition: vmac.h:38
#define begin
Definition: vmac.h:31
HitPattern trackerExpectedHitsOuter_
hit pattern used for expected crossed layers before the first trajectory&#39;s hit
Definition: TrackBase.h:301
void reco::TrackBase::setTrackerExpectedHitsOuter ( const TrackingRecHit hit,
size_t  i 
)
inline

Definition at line 263 of file TrackBase.h.

References reco::HitPattern::set(), and trackerExpectedHitsOuter_.

263 { trackerExpectedHitsOuter_.set( hit, i ); }
int i
Definition: DBlmapReader.cc:9
void set(const I &begin, const I &end)
Definition: HitPattern.h:139
HitPattern trackerExpectedHitsOuter_
hit pattern used for expected crossed layers before the first trajectory&#39;s hit
Definition: TrackBase.h:301
void reco::TrackBase::setTrackerExpectedHitsOuter ( const HitPattern hitP)
inline

Definition at line 268 of file TrackBase.h.

References trackerExpectedHitsOuter_.

HitPattern trackerExpectedHitsOuter_
hit pattern used for expected crossed layers before the first trajectory&#39;s hit
Definition: TrackBase.h:301
double reco::TrackBase::theta ( ) const
inline
double reco::TrackBase::thetaError ( ) const
inline
const HitPattern& reco::TrackBase::trackerExpectedHitsInner ( ) const
inline

Access the hit pattern counting (in the Tracker) the number of expected crossed layers before the first trajectory's hit.

Definition at line 225 of file TrackBase.h.

References trackerExpectedHitsInner_.

Referenced by CheckHitPattern::analyze(), PFCheckHitPattern::analyze(), TkConvValidator::analyze(), IsolatedTracksNxN::analyze(), PrimaryVertexAnalyzer4PU::fillTrackHistos(), spr::goodTrack(), CheckHitPattern::print(), PFCheckHitPattern::print(), IsolatedTracksNxN::printTrack(), and FWConvTrackHitsDetailView::setTextInfo().

225 { return trackerExpectedHitsInner_; }
HitPattern trackerExpectedHitsInner_
hit pattern used for expected crossed layers after the last trajectory&#39;s hit
Definition: TrackBase.h:299
const HitPattern& reco::TrackBase::trackerExpectedHitsOuter ( ) const
inline

Access the hit pattern counting (in the Tracker) the number of expected crossed layers after the last trajectory's hit.

Definition at line 227 of file TrackBase.h.

References trackerExpectedHitsOuter_.

Referenced by IsolatedTracksNxN::analyze(), PrimaryVertexAnalyzer4PU::fillTrackHistos(), spr::goodTrack(), PFDisplacedVertexHelper::isTrackSelected(), and IsolatedTracksNxN::printTrack().

227 { return trackerExpectedHitsOuter_; }
HitPattern trackerExpectedHitsOuter_
hit pattern used for expected crossed layers before the first trajectory&#39;s hit
Definition: TrackBase.h:301
double reco::TrackBase::validFraction ( ) const
inline

fraction of valid hits on the track

Definition at line 236 of file TrackBase.h.

References hitPattern_, reco::HitPattern::numberOfLostTrackerHits(), reco::HitPattern::numberOfValidTrackerHits(), trackerExpectedHitsInner_, trackerExpectedHitsOuter_, and TrackValidation_HighPurity_cff::valid.

236  {
241  if ((valid+lost+lostIn+lostOut)==0) return -1;
242  return valid/(1.0*(valid+lost+lostIn+lostOut));
243  }
int numberOfLostTrackerHits() const
Definition: HitPattern.h:615
HitPattern trackerExpectedHitsInner_
hit pattern used for expected crossed layers after the last trajectory&#39;s hit
Definition: TrackBase.h:299
int numberOfValidTrackerHits() const
Definition: HitPattern.h:558
HitPattern hitPattern_
hit pattern
Definition: TrackBase.h:297
HitPattern trackerExpectedHitsOuter_
hit pattern used for expected crossed layers before the first trajectory&#39;s hit
Definition: TrackBase.h:301
const Point& reco::TrackBase::vertex ( ) const
inline
double reco::TrackBase::vx ( ) const
inline
double reco::TrackBase::vy ( ) const
inline
double reco::TrackBase::vz ( ) const
inline

Member Data Documentation

std::string const TrackBase::algoNames
static
Initial value:
= { "undefAlgorithm", "ctf", "rs", "cosmics",
"iter0", "iter1", "iter2","iter3","iter4","iter5","iter6","iter7","iter8","iter9","iter10",
"outInEcalSeededConv","inOutEcalSeededConv",
"nuclInter",
"standAloneMuon","globalMuon","cosmicStandAloneMuon","cosmicGlobalMuon",
"iter1LargeD0","iter2LargeD0","iter3LargeD0","iter4LargeD0","iter5LargeD0",
"bTagGhostTracks",
"beamhalo" ,
"gsf" }

Definition at line 92 of file TrackBase.h.

Referenced by algoByName(), and algoName().

uint8_t reco::TrackBase::algorithm_
private

track algorithm

Definition at line 318 of file TrackBase.h.

Referenced by algo(), algoName(), and setAlgorithm().

char reco::TrackBase::charge_
private

electric charge

Definition at line 316 of file TrackBase.h.

Referenced by charge().

float reco::TrackBase::chi2_
private

chi-squared

Definition at line 306 of file TrackBase.h.

Referenced by chi2(), and normalizedChi2().

float reco::TrackBase::covariance_[covarianceSize]
private

perigee 5x5 covariance matrix

Definition at line 304 of file TrackBase.h.

Referenced by covariance(), error(), fill(), and TrackBase().

HitPattern reco::TrackBase::hitPattern_
private
Vector reco::TrackBase::momentum_
private

momentum vector at innermost point

Definition at line 310 of file TrackBase.h.

Referenced by eta(), lambda(), momentum(), p(), phi(), pt(), px(), py(), pz(), and theta().

float reco::TrackBase::ndof_
private

number of degrees of freedom

Definition at line 313 of file TrackBase.h.

Referenced by ndof(), and normalizedChi2().

signed char reco::TrackBase::nLoops_
private

number of loops made during the building of the trajectory of a looper particle

Definition at line 323 of file TrackBase.h.

Referenced by isLooper(), nLoops(), and setNLoops().

uint8_t reco::TrackBase::quality_
private

track quality

Definition at line 320 of file TrackBase.h.

Referenced by quality(), qualityMask(), setQuality(), and setQualityMask().

std::string const TrackBase::qualityNames = { "loose", "tight", "highPurity", "confirmed", "goodIterative", "looseSetWithPV", "highPuritySetWithPV"}
static

Definition at line 96 of file TrackBase.h.

Referenced by qualityByName(), and qualityName().

HitPattern reco::TrackBase::trackerExpectedHitsInner_
private

hit pattern used for expected crossed layers after the last trajectory's hit

Definition at line 299 of file TrackBase.h.

Referenced by setTrackerExpectedHitsInner(), trackerExpectedHitsInner(), and validFraction().

HitPattern reco::TrackBase::trackerExpectedHitsOuter_
private

hit pattern used for expected crossed layers before the first trajectory's hit

Definition at line 301 of file TrackBase.h.

Referenced by setTrackerExpectedHitsOuter(), trackerExpectedHitsOuter(), and validFraction().

Point reco::TrackBase::vertex_
private

innermost (reference) point on track

Definition at line 308 of file TrackBase.h.

Referenced by referencePoint(), vertex(), vx(), vy(), and vz().