#include <DataFormats/TrackReco/interface/Track.h>
Public Member Functions | |
const TrackExtraRef & | extra () const |
reference to "extra" object More... | |
CovarianceMatrix & | fillInner (CovarianceMatrix &v) const |
CovarianceMatrix & | fillOuter (CovarianceMatrix &v) const |
fill outermost trajectory state curvilinear errors More... | |
unsigned short | found () const |
Number of valid hits on track. More... | |
unsigned int | innerDetId () const |
DetId of the detector on which surface the innermost state is located. More... | |
const math::XYZVector & | innerMomentum () const |
momentum vector at the innermost hit position More... | |
bool | innerOk () const |
return true if the innermost hit is valid More... | |
const math::XYZPoint & | innerPosition () const |
position of the innermost hit More... | |
CovarianceMatrix | innerStateCovariance () const |
innermost trajectory state curvilinear errors More... | |
unsigned short | lost () const |
Number of lost (=invalid) hits on track. More... | |
unsigned int | outerDetId () const |
DetId of the detector on which surface the outermost state is located. More... | |
double | outerEta () const |
pseudorapidity of the momentum vector at the outermost hit position More... | |
const math::XYZVector & | outerMomentum () const |
momentum vector at the outermost hit position More... | |
bool | outerOk () const |
return true if the outermost hit is valid More... | |
double | outerP () const |
magnitude of momentum vector at the outermost hit position More... | |
double | outerPhi () const |
azimuthal angle of the momentum vector at the outermost hit position More... | |
const math::XYZPoint & | outerPosition () const |
position of the outermost hit More... | |
double | outerPt () const |
transverse momentum at the outermost hit position More... | |
double | outerPx () const |
x coordinate of momentum vector at the outermost hit position More... | |
double | outerPy () const |
y coordinate of momentum vector at the outermost hit position More... | |
double | outerPz () const |
z coordinate of momentum vector at the outermost hit position More... | |
double | outerRadius () const |
polar radius of the outermost hit position More... | |
CovarianceMatrix | outerStateCovariance () const |
outermost trajectory state curvilinear errors More... | |
double | outerTheta () const |
polar angle of the momentum vector at the outermost hit position More... | |
double | outerX () const |
x coordinate of the outermost hit position More... | |
double | outerY () const |
y coordinate of the outermost hit position More... | |
double | outerZ () const |
z coordinate of the outermost hit position More... | |
TrackingRecHitRef | recHit (size_t i) const |
Get i-th hit on the track. More... | |
auto | recHits () const |
Access to reconstructed hits on the track. More... | |
trackingRecHit_iterator | recHitsBegin () const |
Iterator to first hit on the track. More... | |
trackingRecHit_iterator | recHitsEnd () const |
Iterator to last hit on the track. More... | |
size_t | recHitsSize () const |
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits). More... | |
const TrackResiduals & | residuals () const |
get the residuals More... | |
const PropagationDirection & | seedDirection () const |
direction of how the hits were sorted in the original seed More... | |
const edm::RefToBase< TrajectorySeed > & | seedRef () const |
void | setExtra (const TrackExtraRef &ref) |
set reference to "extra" object More... | |
Track () | |
default constructor More... | |
Track (double chi2, double ndof, const Point &referencePoint, const Vector &momentum, int charge, const CovarianceMatrix &, TrackAlgorithm=undefAlgorithm, TrackQuality quality=undefQuality, float t0=0, float beta=0, float covt0t0=-1., float covbetabeta=-1.) | |
constructor from fit parameters and error matrix More... | |
~Track () override | |
virtual destructor More... | |
Public Member Functions inherited from reco::TrackBase | |
TrackAlgorithm | algo () const |
AlgoMask | algoMask () const |
unsigned long long | algoMaskUL () const |
std::string | algoName () const |
bool | appendHitPattern (const TrackingRecHit &hit, const TrackerTopology &ttopo) |
append a single hit to the HitPattern More... | |
bool | appendHitPattern (const DetId &id, TrackingRecHit::Type hitType, const TrackerTopology &ttopo) |
bool | appendHitPattern (const uint16_t pattern, TrackingRecHit::Type hitType) |
template<typename C > | |
bool | appendHits (const C &c, const TrackerTopology &ttopo) |
append hit patterns from vector of hit references More... | |
template<typename I > | |
bool | appendHits (const I &begin, const I &end, const TrackerTopology &ttopo) |
bool | appendMuonHitPattern (const DetId &id, TrackingRecHit::Type hitType) |
bool | appendTrackerHitPattern (uint16_t subdet, uint16_t layer, uint16_t stereo, TrackingRecHit::Type hitType) |
double | beta () const |
velocity at the reference point in natural units More... | |
double | betaError () const |
error on beta More... | |
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 covariance matrix (i, j = 0, ... 4) More... | |
double | covBetaBeta () const |
error on beta More... | |
double | covt0t0 () const |
error on t0 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). 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 | dxyError (Point const &vtx, math::Error< 3 >::type const &vertexCov) const |
error on dxy with respect to a user-given reference point + uncertainty (i.e. reco::Vertex position) More... | |
double | dxyError (const BeamSpot &theBeamSpot) const |
error on dxy with respect to a user-given beamspot 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) 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... | |
CovarianceMatrix & | fill (CovarianceMatrix &v) const |
fill SMatrix More... | |
const HitPattern & | hitPattern () const |
Access the hit pattern, indicating in which Tracker layers the track has hits. More... | |
bool | isAlgoInMask (TrackAlgorithm a) const |
bool | isLooper () const |
bool | isTimeOk () const |
return true if timing measurement is usable More... | |
double | lambda () const |
Lambda angle. More... | |
double | lambdaError () const |
error on lambda More... | |
int | missingInnerHits () const |
number of hits expected from inner track extrapolation but missing More... | |
int | missingOuterHits () const |
number of hits expected from outer track extrapolation but missing More... | |
const Vector & | momentum () 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... | |
TrackAlgorithm | originalAlgo () const |
double | p () const |
momentum vector magnitude More... | |
double | p2 () const |
momentum vector magnitude square 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 | pt2 () const |
track transverse momentum square More... | |
double | ptError () const |
error on Pt (set to 1000 TeV if charge==0 for safety) More... | |
double | ptError2 () const |
error on Pt (set to 1000**2 TeV**2 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 Point & | referencePoint () const |
Reference point on the track. More... | |
void | resetHitPattern () |
Sets HitPattern as empty. More... | |
void | setAlgoMask (AlgoMask a) |
void | setAlgorithm (const TrackAlgorithm a) |
Track algorithm. More... | |
void | setNLoops (signed char value) |
void | setOriginalAlgorithm (const TrackAlgorithm a) |
void | setQuality (const TrackQuality) |
void | setQualityMask (int qualMask) |
void | setStopReason (uint8_t value) |
uint8_t | stopReason () const |
double | t0 () const |
time at the reference point More... | |
double | t0Error () const |
error on t0 More... | |
double | theta () const |
polar angle More... | |
double | thetaError () const |
error on theta More... | |
TrackBase () | |
default constructor More... | |
TrackBase (double chi2, double ndof, const Point &vertex, const Vector &momentum, int charge, const CovarianceMatrix &cov, TrackAlgorithm=undefAlgorithm, TrackQuality quality=undefQuality, signed char nloops=0, uint8_t stopReason=0, float t0=0.f, float beta=0.f, float covt0t0=-1.f, float covbetabeta=-1.f) | |
constructor from fit parameters and error matrix More... | |
double | validFraction () const |
fraction of valid hits on the track More... | |
const Point & | vertex () 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... | |
virtual | ~TrackBase () |
virtual destructor More... | |
Private Attributes | |
TrackExtraRef | extra_ |
Reference to additional information stored only on RECO. More... | |
This class describes the reconstructed tracks that are stored in the AOD and RECO. It also contains a reference to more detailed information about each track, that is stoed in the TrackExtra object, available only in RECO.
Note that most of the functions provided in this Track class rely on the existance of the TrackExtra object, so will not work on AOD.
The most useful functions are those provided in the TrackBase class from which this inherits, all of which work on AOD.
Track::Track | ( | double | chi2, |
double | ndof, | ||
const Point & | referencePoint, | ||
const Vector & | momentum, | ||
int | charge, | ||
const CovarianceMatrix & | cov, | ||
TrackAlgorithm | algo = undefAlgorithm , |
||
TrackQuality | quality = undefQuality , |
||
float | t0 = 0 , |
||
float | beta = 0 , |
||
float | covt0t0 = -1. , |
||
float | covbetabeta = -1. |
||
) |
constructor from fit parameters and error matrix
Definition at line 5 of file Track.cc.
|
inline |
reference to "extra" object
Definition at line 139 of file Track.h.
References extra_.
Referenced by BeamHaloAnalyzer::analyze(), PFElecTkProducer::applySelection(), FWTrackProxyBuilder::build(), FWConvTrackHitsDetailView::build(), ConversionProducer::buildCollection(), CSCHaloAlgo::Calculate(), ConversionProducer::checkPhi(), helper::GsfElectronCollectionStoreManager::cloneAndStore(), PFElecTkProducer::findPfRef(), AlCaHOCalibProducer::getFreeTrajectoryState(), FWPFTrackUtils::getTrack(), trajectoryStateTransform::innerFreeState(), trajectoryStateTransform::innerStateOnSurface(), MultiTrajectoryStateTransform::innerStateOnSurface(), trajectoryStateTransform::outerFreeState(), trajectoryStateTransform::outerStateOnSurface(), MultiTrajectoryStateTransform::outerStateOnSurface(), MultiTrackSelector::processMVA(), CosmicTrackSelector::produce(), reco::modules::TrackFullCloneSelectorBase< Selector >::produce(), TrackExtrapolator::propagateTrackToVolume(), AnalyticalTrackSelector::run(), and MultiTrackSelector::select().
|
inline |
Definition at line 76 of file Track.h.
References extra_, and findQualityFiles::v.
|
inline |
fill outermost trajectory state curvilinear errors
Definition at line 74 of file Track.h.
References extra_, and findQualityFiles::v.
|
inline |
Number of valid hits on track.
Definition at line 142 of file Track.h.
References reco::TrackBase::numberOfValidHits().
Referenced by lowptgsfeleid::features_V1(), ConversionProducer::trackQualityFilter(), and TSGForOIDNN::updateFeatureMap().
|
inline |
DetId of the detector on which surface the innermost state is located.
Definition at line 82 of file Track.h.
References extra_.
Referenced by DTChamberEfficiency::analyze(), helper::GsfElectronCollectionStoreManager::cloneAndStore(), CalibrationTrackSelectorFromDetIdList::makeCandidate(), reco::modules::CosmicTrackSplitter::makeCandidate(), reco::modules::TrackerTrackHitFilter::makeCandidate(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), TrackListCombiner::produce(), MuonTrackProducer::produce(), L3TkMuonProducer::produce(), CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::TrackFullCloneSelectorBase< Selector >::produce(), AnalyticalTrackSelector::run(), and GlobalMuonRefitter::transform().
|
inline |
momentum vector at the innermost hit position
Definition at line 59 of file Track.h.
References extra_.
Referenced by ConversionProducer::buildCollection(), ConversionProducer::checkPhi(), helper::GsfElectronCollectionStoreManager::cloneAndStore(), PFEGammaProducer::createSingleLegConversions(), AlCaHOCalibProducer::getFreeTrajectoryState(), trajectoryStateTransform::innerFreeState(), ConversionProducer::preselectTrackPair(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), TrackListCombiner::produce(), MuonTrackProducer::produce(), L3TkMuonProducer::produce(), CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::TrackFullCloneSelectorBase< Selector >::produce(), spr::propagateCosmicCALO(), MuonTrackingRegionByPtBuilder::region(), MuonTrackingRegionBuilder::region(), and AnalyticalTrackSelector::run().
|
inline |
return true if the innermost hit is valid
Definition at line 53 of file Track.h.
References extra_, edm::Ref< C, T, F >::isAvailable(), and edm::Ref< C, T, F >::isNonnull().
Referenced by helper::GsfElectronCollectionStoreManager::cloneAndStore(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), TrackListCombiner::produce(), MuonTrackProducer::produce(), L3TkMuonProducer::produce(), CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::TrackFullCloneSelectorBase< Selector >::produce(), and AnalyticalTrackSelector::run().
|
inline |
position of the innermost hit
Definition at line 56 of file Track.h.
References extra_.
Referenced by BeamHaloAnalyzer::analyze(), ConversionProducer::buildCollection(), TrackEfficiencyMonitor::checkSemiCylinder(), helper::GsfElectronCollectionStoreManager::cloneAndStore(), PFEGammaProducer::createSingleLegConversions(), AlignmentTrackSelector::detailedHitsCheck(), AlCaHOCalibProducer::getFreeTrajectoryState(), trajectoryStateTransform::innerFreeState(), CalibrationTrackSelectorFromDetIdList::makeCandidate(), reco::modules::CosmicTrackSplitter::makeCandidate(), reco::modules::TrackerTrackHitFilter::makeCandidate(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), TrackListCombiner::produce(), MuonTrackProducer::produce(), L3TkMuonProducer::produce(), CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::TrackFullCloneSelectorBase< Selector >::produce(), spr::propagateCosmicCALO(), AnalyticalTrackSelector::run(), PropagateToMuon::startingState(), and SiStripFineDelayTOF::trackParameters().
|
inline |
innermost trajectory state curvilinear errors
Definition at line 71 of file Track.h.
References extra_.
Referenced by helper::GsfElectronCollectionStoreManager::cloneAndStore(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), TrackListCombiner::produce(), MuonTrackProducer::produce(), L3TkMuonProducer::produce(), CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::TrackFullCloneSelectorBase< Selector >::produce(), and AnalyticalTrackSelector::run().
|
inline |
Number of lost (=invalid) hits on track.
Definition at line 145 of file Track.h.
References reco::TrackBase::numberOfLostHits().
|
inline |
DetId of the detector on which surface the outermost state is located.
Definition at line 79 of file Track.h.
References extra_.
Referenced by helper::GsfElectronCollectionStoreManager::cloneAndStore(), CalibrationTrackSelectorFromDetIdList::makeCandidate(), reco::modules::CosmicTrackSplitter::makeCandidate(), reco::modules::TrackerTrackHitFilter::makeCandidate(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), TrackListCombiner::produce(), MuonTrackProducer::produce(), L3TkMuonProducer::produce(), CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::TrackFullCloneSelectorBase< Selector >::produce(), AnalyticalTrackSelector::run(), and GlobalMuonRefitter::transform().
|
inline |
pseudorapidity of the momentum vector at the outermost hit position
Definition at line 127 of file Track.h.
References extra_.
Referenced by EwkMuLumiMonitorDQM::tkIso().
|
inline |
momentum vector at the outermost hit position
Definition at line 65 of file Track.h.
References extra_.
Referenced by ConversionProducer::buildCollection(), CSCHaloAlgo::Calculate(), helper::GsfElectronCollectionStoreManager::cloneAndStore(), PFEGammaProducer::createSingleLegConversions(), trajectoryStateTransform::outerFreeState(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), TrackListCombiner::produce(), MuonTrackProducer::produce(), L3TkMuonProducer::produce(), CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::TrackFullCloneSelectorBase< Selector >::produce(), spr::propagateCosmicCALO(), AnalyticalTrackSelector::run(), and SiStripFineDelayTOF::trackParameters().
|
inline |
return true if the outermost hit is valid
Definition at line 50 of file Track.h.
References extra_, edm::Ref< C, T, F >::isAvailable(), and edm::Ref< C, T, F >::isNonnull().
Referenced by helper::GsfElectronCollectionStoreManager::cloneAndStore(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), TrackListCombiner::produce(), MuonTrackProducer::produce(), L3TkMuonProducer::produce(), CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::TrackFullCloneSelectorBase< Selector >::produce(), and AnalyticalTrackSelector::run().
|
inline |
magnitude of momentum vector at the outermost hit position
Definition at line 118 of file Track.h.
References extra_.
|
inline |
azimuthal angle of the momentum vector at the outermost hit position
Definition at line 124 of file Track.h.
References extra_.
Referenced by EwkMuLumiMonitorDQM::tkIso().
|
inline |
position of the outermost hit
Definition at line 62 of file Track.h.
References extra_.
Referenced by IsolatedTracksNxN::analyze(), BeamHaloAnalyzer::analyze(), helper::GsfElectronCollectionStoreManager::cloneAndStore(), PFEGammaProducer::createSingleLegConversions(), AlignmentTrackSelector::detailedHitsCheck(), CSCEfficiency::filter(), CalibrationTrackSelectorFromDetIdList::makeCandidate(), reco::modules::CosmicTrackSplitter::makeCandidate(), reco::modules::TrackerTrackHitFilter::makeCandidate(), trajectoryStateTransform::outerFreeState(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), TrackListCombiner::produce(), MuonTrackProducer::produce(), L3TkMuonProducer::produce(), CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::TrackFullCloneSelectorBase< Selector >::produce(), spr::propagateCosmicCALO(), AnalyticalTrackSelector::run(), PropagateToMuon::startingState(), and SiStripFineDelayTOF::trackParameters().
|
inline |
transverse momentum at the outermost hit position
Definition at line 121 of file Track.h.
References extra_.
|
inline |
x coordinate of momentum vector at the outermost hit position
Definition at line 100 of file Track.h.
References extra_.
Referenced by AlCaHOCalibProducer::getFreeTrajectoryState(), and TrackExtrapolator::propagateTrackToVolume().
|
inline |
y coordinate of momentum vector at the outermost hit position
Definition at line 103 of file Track.h.
References extra_.
Referenced by AlCaHOCalibProducer::getFreeTrajectoryState(), and TrackExtrapolator::propagateTrackToVolume().
|
inline |
z coordinate of momentum vector at the outermost hit position
Definition at line 106 of file Track.h.
References extra_.
Referenced by AlCaHOCalibProducer::getFreeTrajectoryState(), and TrackExtrapolator::propagateTrackToVolume().
|
inline |
polar radius of the outermost hit position
Definition at line 133 of file Track.h.
References extra_.
|
inline |
outermost trajectory state curvilinear errors
Definition at line 68 of file Track.h.
References extra_.
Referenced by helper::GsfElectronCollectionStoreManager::cloneAndStore(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), TrackListCombiner::produce(), MuonTrackProducer::produce(), L3TkMuonProducer::produce(), CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::TrackFullCloneSelectorBase< Selector >::produce(), and AnalyticalTrackSelector::run().
|
inline |
polar angle of the momentum vector at the outermost hit position
Definition at line 130 of file Track.h.
References extra_.
|
inline |
x coordinate of the outermost hit position
Definition at line 109 of file Track.h.
References extra_.
Referenced by AlCaHOCalibProducer::getFreeTrajectoryState(), and TrackExtrapolator::propagateTrackToVolume().
|
inline |
y coordinate of the outermost hit position
Definition at line 112 of file Track.h.
References extra_.
Referenced by AlCaHOCalibProducer::getFreeTrajectoryState(), and TrackExtrapolator::propagateTrackToVolume().
|
inline |
z coordinate of the outermost hit position
Definition at line 115 of file Track.h.
References extra_.
Referenced by AlCaHOCalibProducer::getFreeTrajectoryState(), and TrackExtrapolator::propagateTrackToVolume().
|
inline |
Get i-th hit on the track.
Definition at line 94 of file Track.h.
References extra_, and mps_fire::i.
Referenced by OutsideInMuonSeeder::doDebug(), MuonAlignmentAnalyzer::doMatching(), SeedGeneratorFromProtoTracksEDProducer::produce(), TrackListCombiner::produce(), reco::TransientTrack::recHit(), and SeedFromProtoTrack::SeedFromProtoTrack().
|
inline |
Access to reconstructed hits on the track.
Definition at line 85 of file Track.h.
References extra_.
Referenced by TrackSplittingMonitor::analyze(), CosmicSplitterValidation::analyze(), PrimaryVertexValidation::analyze(), SegmentsTrackAssociator::associate(), helper::GsfElectronCollectionStoreManager::cloneAndStore(), CalibrationTrackSelector::detailedHitsCheck(), AlignmentTrackSelector::detailedHitsCheck(), ConvBremSeedProducer::isGsfTrack(), spr::matchedSimTrack(), MuonResidualsFromTrack::MuonResidualsFromTrack(), and reco::modules::TrackFullCloneSelectorBase< Selector >::produce().
|
inline |
Iterator to first hit on the track.
Definition at line 88 of file Track.h.
References extra_.
Referenced by RPCRecHitProbability::analyze(), RPCMonitorDigi::analyze(), TestOutliers::analyze(), FWTrackTrackingRecHitProxyBuilder::build(), SiStripFineDelayHit::closestCluster(), TrackQuality::evaluate(), ConversionHitChecker::nSharedHits(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), MuonTrackProducer::produce(), L3TkMuonProducer::produce(), CosmicTrackSelector::produce(), TrackerCleaner< T >::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::CosmicTrackSplitter::produce(), reco::modules::TrackerTrackHitFilter::produceFromTrack(), SiStripFineDelayHit::rechit(), reco::TransientTrack::recHitsBegin(), AnalyticalTrackSelector::run(), CosmicTrackSelector::select(), and CosmicMuonLinksProducer::sharedHits().
|
inline |
Iterator to last hit on the track.
Definition at line 91 of file Track.h.
References extra_.
Referenced by RPCRecHitProbability::analyze(), RPCMonitorDigi::analyze(), TestOutliers::analyze(), FWTrackTrackingRecHitProxyBuilder::build(), SiStripFineDelayHit::closestCluster(), SiStripFineDelayHit::detId(), TrackQuality::evaluate(), ConversionHitChecker::nSharedHits(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), MuonTrackProducer::produce(), L3TkMuonProducer::produce(), CosmicTrackSelector::produce(), TrackerCleaner< T >::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::CosmicTrackSplitter::produce(), reco::modules::TrackerTrackHitFilter::produceFromTrack(), SiStripFineDelayHit::rechit(), reco::TransientTrack::recHitsEnd(), AnalyticalTrackSelector::run(), CosmicTrackSelector::select(), and CosmicMuonLinksProducer::sharedHits().
|
inline |
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits).
Definition at line 97 of file Track.h.
References extra_.
Referenced by helper::GsfElectronCollectionStoreManager::cloneAndStore(), OutsideInMuonSeeder::doDebug(), MultiTrackSelector::processMVA(), SeedGeneratorFromProtoTracksEDProducer::produce(), TrackListCombiner::produce(), reco::TransientTrack::recHitsSize(), SeedFromProtoTrack::SeedFromProtoTrack(), and MultiTrackSelector::select().
|
inline |
get the residuals
Definition at line 158 of file Track.h.
References extra_.
Referenced by CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), and AnalyticalTrackSelector::run().
|
inline |
direction of how the hits were sorted in the original seed
Definition at line 148 of file Track.h.
References extra_.
Referenced by helper::GsfElectronCollectionStoreManager::cloneAndStore(), CalibrationTrackSelectorFromDetIdList::makeCandidate(), reco::modules::CosmicTrackSplitter::makeCandidate(), reco::modules::TrackerTrackHitFilter::makeCandidate(), helper::MuonCollectionStoreManager::processMuon(), helper::TrackCollectionStoreManager::processTrack(), TrackListCombiner::produce(), MuonTrackProducer::produce(), CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), reco::modules::TrackFullCloneSelectorBase< Selector >::produce(), and AnalyticalTrackSelector::run().
|
inline |
return the edm::reference to the trajectory seed in the original seeds collection. If the collection has been dropped from the Event, the reference may be invalid. Its validity should be tested, before the reference is actually used.
Definition at line 155 of file Track.h.
References extra_.
Referenced by TestOutliers::analyze(), PFElecTkProducer::applySelection(), PFElecTkProducer::findPfRef(), CalibrationTrackSelectorFromDetIdList::makeCandidate(), reco::modules::CosmicTrackSplitter::makeCandidate(), reco::modules::TrackerTrackHitFilter::makeCandidate(), TrackListCombiner::produce(), CosmicTrackSelector::produce(), reco::modules::HICaloCompatibleTrackSelector::produce(), and AnalyticalTrackSelector::run().
|
inline |
set reference to "extra" object
Definition at line 136 of file Track.h.
References extra_.
Referenced by helper::GsfElectronCollectionStoreManager::cloneAndStore(), helper::MuonCollectionStoreManager::processMuon(), SingleLongTrackProducer::produce(), and SeedToTrackProducer::produce().
|
private |
Reference to additional information stored only on RECO.
Definition at line 162 of file Track.h.
Referenced by extra(), fillInner(), fillOuter(), innerDetId(), innerMomentum(), innerOk(), innerPosition(), innerStateCovariance(), outerDetId(), outerEta(), outerMomentum(), outerOk(), outerP(), outerPhi(), outerPosition(), outerPt(), outerPx(), outerPy(), outerPz(), outerRadius(), outerStateCovariance(), outerTheta(), outerX(), outerY(), outerZ(), recHit(), recHits(), recHitsBegin(), recHitsEnd(), recHitsSize(), residuals(), seedDirection(), seedRef(), and setExtra().