CMS 3D CMS Logo

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

#include <TrackTransientTrack.h>

Inheritance diagram for reco::TrackTransientTrack:
reco::Track reco::BasicTransientTrack reco::TrackBase BasicReferenceCounted

Public Member Functions

TrackCharge charge () const
 
const MagneticFieldfield () const
 
TrajectoryStateOnSurface impactPointState () const
 
bool impactPointStateAvailable () const
 
TrajectoryStateClosestToPoint impactPointTSCP () const
 
FreeTrajectoryState initialFreeState () const
 
TrajectoryStateOnSurface innermostMeasurementState () const
 
TrackTransientTrackoperator= (const TrackTransientTrack &tt)
 
TrajectoryStateOnSurface outermostMeasurementState () const
 
TrackRef persistentTrackRef () const
 
void setBeamSpot (const reco::BeamSpot &beamSpot)
 
void setES (const edm::EventSetup &)
 
void setTrackingGeometry (const edm::ESHandle< GlobalTrackingGeometry > &)
 
TrajectoryStateClosestToBeamLine stateAtBeamLine () const
 
TrajectoryStateOnSurface stateOnSurface (const GlobalPoint &point) const
 
const Tracktrack () const
 
TrackBaseRef trackBaseRef () const
 
 TrackTransientTrack ()
 
 TrackTransientTrack (const Track &tk, const MagneticField *field)
 
 TrackTransientTrack (const TrackRef &tk, const MagneticField *field)
 
 TrackTransientTrack (const TrackRef &tk, const MagneticField *field, const edm::ESHandle< GlobalTrackingGeometry > &trackingGeometry)
 
 TrackTransientTrack (const Track &tk, const MagneticField *field, const edm::ESHandle< GlobalTrackingGeometry > &trackingGeometry)
 
 TrackTransientTrack (const TrackTransientTrack &tt)
 
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint (const GlobalPoint &point) const
 
- Public Member Functions inherited from reco::Track
const TrackExtraRefextra () const
 reference to "extra" object More...
 
CovarianceMatrixfillInner (CovarianceMatrix &v) const
 
CovarianceMatrixfillOuter (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::XYZVectorinnerMomentum () const
 momentum vector at the innermost hit position More...
 
bool innerOk () const
 return true if the innermost hit is valid More...
 
const math::XYZPointinnerPosition () 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::XYZVectorouterMomentum () 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::XYZPointouterPosition () 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...
 
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 TrackResidualsresiduals () const
 
double residualX (int position) const
 
double residualY (int position) const
 
PropagationDirection seedDirection () const
 direction of how the hits were sorted in the original seed More...
 
edm::RefToBase< TrajectorySeedseedRef () 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)
 constructor from fit parameters and error matrix More...
 
virtual ~Track ()
 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)
 
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)
 
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 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 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...
 
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 isAlgoInMask (TrackAlgorithm a) const
 
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...
 
TrackAlgorithm originalAlgo () const
 
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 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)
 
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)
 constructor from fit parameters and error matrix 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...
 
virtual ~TrackBase ()
 virtual destructor More...
 
- Public Member Functions inherited from reco::BasicTransientTrack
virtual CandidatePtr candidate () const
 
virtual ~BasicTransientTrack ()
 
- Public Member Functions inherited from BasicReferenceCounted
void addReference () const
 
 BasicReferenceCounted ()
 
 BasicReferenceCounted (const BasicReferenceCounted &)
 
const BasicReferenceCountedoperator= (const BasicReferenceCounted &)
 
unsigned int references () const
 
void removeReference () const
 
virtual ~BasicReferenceCounted ()
 

Private Types

enum  CacheStates { kUnset, kSetting, kSet }
 

Private Attributes

TSCPBuilderNoMaterial builder
 
FreeTrajectoryState initialFTS
 
TrajectoryStateClosestToPoint initialTSCP
 
TrajectoryStateOnSurface initialTSOS
 
std::atomic< char > m_SCTBL
 
std::atomic< char > m_TSCP
 
std::atomic< char > m_TSOS
 
reco::BeamSpot theBeamSpot
 
const MagneticFieldtheField
 
edm::ESHandle
< GlobalTrackingGeometry
theTrackingGeometry
 
TrackRef tkr_
 
TrajectoryStateClosestToBeamLine trajectoryStateClosestToBeamLine
 

Additional Inherited Members

- Public Types inherited from reco::TrackBase
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 std::bitset< algoSizeAlgoMask
 algo mask 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,
  initialStep = 4, lowPtTripletStep = 5, pixelPairStep = 6, detachedTripletStep = 7,
  mixedTripletStep = 8, pixelLessStep = 9, tobTecStep = 10, jetCoreRegionalStep = 11,
  conversionStep = 12, muonSeededStepInOut = 13, muonSeededStepOutIn = 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, hltPixel = 30, hltIter0 = 31,
  hltIter1 = 32, hltIter2 = 33, hltIter3 = 34, hltIter4 = 35,
  hltIterX = 36, hiRegitMuInitialStep = 37, hiRegitMuLowPtTripletStep = 38, hiRegitMuPixelPairStep = 39,
  hiRegitMuDetachedTripletStep = 40, hiRegitMuMixedTripletStep = 41, hiRegitMuPixelLessStep = 42, hiRegitMuTobTecStep = 43,
  hiRegitMuMuonSeededStepInOut = 44, hiRegitMuMuonSeededStepOutIn = 45, algoSize = 46
}
 track algorithm More...
 
enum  TrackQuality {
  undefQuality = -1, loose = 0, tight = 1, highPurity = 2,
  confirmed = 3, goodIterative = 4, looseSetWithPV = 5, highPuritySetWithPV = 6,
  discarded = 7, qualitySize = 8
}
 track quality More...
 
typedef math::XYZVector Vector
 spatial vector More...
 
- Public Types inherited from reco::BasicTransientTrack
typedef BasicTransientTrack BTT
 
typedef ProxyBase< BTT,
CopyUsingClone< BTT > > 
Proxy
 
typedef
ReferenceCountingPointer
< BasicTransientTrack
RCPtr
 
- Static Public Member Functions inherited from reco::TrackBase
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 inherited from reco::TrackBase
static const std::string algoNames []
 
static const std::string qualityNames []
 

Detailed Description

Definition at line 18 of file TrackTransientTrack.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

TrackTransientTrack::TrackTransientTrack ( )

Definition at line 23 of file TrackTransientTrack.cc.

23  :
25 {
26 }
std::atomic< char > m_TSOS
Track()
default constructor
Definition: Track.h:34
const MagneticField * theField
std::atomic< char > m_TSCP
std::atomic< char > m_SCTBL
TrackTransientTrack::TrackTransientTrack ( const Track tk,
const MagneticField field 
)

Definition at line 28 of file TrackTransientTrack.cc.

References trajectoryStateTransform::initialFreeState(), and initialFTS.

28  :
30 {
31 
33 }
std::atomic< char > m_TSOS
Track()
default constructor
Definition: Track.h:34
const MagneticField * theField
std::atomic< char > m_TSCP
std::atomic< char > m_SCTBL
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
FreeTrajectoryState initialFTS
TrackTransientTrack::TrackTransientTrack ( const TrackRef tk,
const MagneticField field 
)

Definition at line 36 of file TrackTransientTrack.cc.

References trajectoryStateTransform::initialFreeState(), and initialFTS.

36  :
37  Track(*tk), tkr_(tk), theField(field), m_TSOS(kUnset), m_TSCP(kUnset), m_SCTBL(kUnset)
38 {
39 
41 }
std::atomic< char > m_TSOS
Track()
default constructor
Definition: Track.h:34
const MagneticField * theField
std::atomic< char > m_TSCP
std::atomic< char > m_SCTBL
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
FreeTrajectoryState initialFTS
TrackTransientTrack::TrackTransientTrack ( const TrackRef tk,
const MagneticField field,
const edm::ESHandle< GlobalTrackingGeometry > &  trackingGeometry 
)

Definition at line 50 of file TrackTransientTrack.cc.

References trajectoryStateTransform::initialFreeState(), and initialFTS.

50  :
52 {
53 
55 }
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
std::atomic< char > m_TSOS
Track()
default constructor
Definition: Track.h:34
const MagneticField * theField
std::atomic< char > m_TSCP
std::atomic< char > m_SCTBL
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
FreeTrajectoryState initialFTS
TrackTransientTrack::TrackTransientTrack ( const Track tk,
const MagneticField field,
const edm::ESHandle< GlobalTrackingGeometry > &  trackingGeometry 
)

Definition at line 43 of file TrackTransientTrack.cc.

References trajectoryStateTransform::initialFreeState(), and initialFTS.

43  :
45 {
46 
48 }
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
std::atomic< char > m_TSOS
Track()
default constructor
Definition: Track.h:34
const MagneticField * theField
std::atomic< char > m_TSCP
std::atomic< char > m_SCTBL
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
FreeTrajectoryState initialFTS
TrackTransientTrack::TrackTransientTrack ( const TrackTransientTrack tt)

Definition at line 58 of file TrackTransientTrack.cc.

References impactPointState(), impactPointTSCP(), initialTSCP, initialTSOS, kSet, m_TSCP, and m_TSOS.

58  :
59  Track(tt), tkr_(tt.persistentTrackRef()), theField(tt.field()),
61 {
62  // see ThreadSafe statement above about the order of operator= and store
63  if (kSet == tt.m_TSOS.load()) {
65  m_TSOS.store(kSet);
66  }
67  // see ThreadSafe statement above about the order of operator= and store
68  if (kSet == tt.m_TSCP.load()) {
70  m_TSCP.store(kSet);
71  }
72 }
const MagneticField * field() const
TrajectoryStateOnSurface impactPointState() const
TrajectoryStateClosestToPoint impactPointTSCP() const
TrajectoryStateOnSurface initialTSOS
std::atomic< char > m_TSOS
Track()
default constructor
Definition: Track.h:34
TrackRef persistentTrackRef() const
FreeTrajectoryState initialFreeState() const
const MagneticField * theField
std::atomic< char > m_TSCP
TrajectoryStateClosestToPoint initialTSCP
FreeTrajectoryState initialFTS

Member Function Documentation

TrackCharge reco::TrackTransientTrack::charge ( void  ) const
inlinevirtual

Implements reco::BasicTransientTrack.

Definition at line 68 of file TrackTransientTrack.h.

References reco::TrackBase::charge().

68 {return Track::charge();}
int charge() const
track electric charge
Definition: TrackBase.h:554
const MagneticField* reco::TrackTransientTrack::field ( ) const
inlinevirtual

Implements reco::BasicTransientTrack.

Definition at line 70 of file TrackTransientTrack.h.

References theField.

70 {return theField;}
const MagneticField * theField
TrajectoryStateOnSurface TrackTransientTrack::impactPointState ( ) const
virtual

Implements reco::BasicTransientTrack.

Definition at line 92 of file TrackTransientTrack.cc.

References TransverseImpactPointExtrapolator::extrapolate(), initialFTS, initialTSOS, kSet, kSetting, kUnset, m_TSOS, FreeTrajectoryState::position(), theField, and tmp.

Referenced by TrackTransientTrack().

93 {
94  // see ThreadSafe statement above about the order of operator= and store
95  if(kSet == m_TSOS.load()) return initialTSOS;
97  auto tmp = tipe.extrapolate(initialFTS, initialFTS.position());
98  char expected = kUnset;
99  if(m_TSOS.compare_exchange_strong(expected, kSetting)) {
100  initialTSOS = tmp;
101  m_TSOS.store(kSet);
102  return initialTSOS;
103  }
104  return tmp;
105 }
TrajectoryStateOnSurface initialTSOS
std::atomic< char > m_TSOS
const MagneticField * theField
GlobalPoint position() const
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
FreeTrajectoryState initialFTS
bool reco::TrackTransientTrack::impactPointStateAvailable ( ) const
inlinevirtual

Implements reco::BasicTransientTrack.

Definition at line 59 of file TrackTransientTrack.h.

References kSet, and m_TSOS.

59 {return (m_TSOS.load()==kSet ? true : false); }
std::atomic< char > m_TSOS
TrajectoryStateClosestToPoint TrackTransientTrack::impactPointTSCP ( ) const
virtual

Implements reco::BasicTransientTrack.

Definition at line 107 of file TrackTransientTrack.cc.

References builder, initialFTS, initialTSCP, kSet, kSetting, kUnset, m_TSCP, FreeTrajectoryState::position(), and tmp.

Referenced by TrackTransientTrack().

108 {
109  // see ThreadSafe statement above about the order of operator= and store
110  if(kSet == m_TSCP.load()) return initialTSCP;
112  char expected = kUnset;
113  if(m_TSCP.compare_exchange_strong(expected, kSetting)) {
114  initialTSCP = tmp;
115  m_TSCP.store(kSet);
116  return initialTSCP;
117  }
118  return tmp;
119 }
TSCPBuilderNoMaterial builder
GlobalPoint position() const
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
std::atomic< char > m_TSCP
TrajectoryStateClosestToPoint initialTSCP
FreeTrajectoryState initialFTS
FreeTrajectoryState reco::TrackTransientTrack::initialFreeState ( ) const
inlinevirtual

Implements reco::BasicTransientTrack.

Definition at line 40 of file TrackTransientTrack.h.

References initialFTS.

40 {return initialFTS;}
FreeTrajectoryState initialFTS
TrajectoryStateOnSurface TrackTransientTrack::innermostMeasurementState ( ) const
virtual

Implements reco::BasicTransientTrack.

Definition at line 127 of file TrackTransientTrack.cc.

References trajectoryStateTransform::innerStateOnSurface(), theField, and theTrackingGeometry.

128 {
129 
131 }
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
const MagneticField * theField
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
TrackTransientTrack& reco::TrackTransientTrack::operator= ( const TrackTransientTrack tt)
TrajectoryStateOnSurface TrackTransientTrack::outermostMeasurementState ( ) const
virtual

Implements reco::BasicTransientTrack.

Definition at line 121 of file TrackTransientTrack.cc.

References trajectoryStateTransform::outerStateOnSurface(), theField, and theTrackingGeometry.

122 {
123 
125 }
TrajectoryStateOnSurface outerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
const MagneticField * theField
TrackRef reco::TrackTransientTrack::persistentTrackRef ( ) const
inline

access to original persistent track

Definition at line 64 of file TrackTransientTrack.h.

References tkr_.

Referenced by ConvertedPhotonProducer::buildCollections(), KinematicVertex::operator reco::Vertex(), and ConversionTrackPairFinder::run().

64 { return tkr_; }
void TrackTransientTrack::setBeamSpot ( const reco::BeamSpot beamSpot)
virtual
void TrackTransientTrack::setES ( const edm::EventSetup setup)
virtual

Implements reco::BasicTransientTrack.

Definition at line 74 of file TrackTransientTrack.cc.

References edm::EventSetup::get(), and theTrackingGeometry.

74  {
75 
77 
78 }
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
const T & get() const
Definition: EventSetup.h:56
void TrackTransientTrack::setTrackingGeometry ( const edm::ESHandle< GlobalTrackingGeometry > &  tg)
virtual

Implements reco::BasicTransientTrack.

Definition at line 80 of file TrackTransientTrack.cc.

References theTrackingGeometry.

80  {
81 
83 
84 }
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
TrajectoryStateClosestToBeamLine TrackTransientTrack::stateAtBeamLine ( ) const
virtual

Implements reco::BasicTransientTrack.

Definition at line 140 of file TrackTransientTrack.cc.

References initialFTS, kSet, kSetting, kUnset, m_SCTBL, theBeamSpot, tmp, and trajectoryStateClosestToBeamLine.

141 {
142  // see ThreadSafe statement above about the order of operator= and store
143  if(kSet == m_SCTBL.load()) return trajectoryStateClosestToBeamLine;
144  TSCBLBuilderNoMaterial blsBuilder;
145  const auto tmp = blsBuilder(initialFTS, theBeamSpot);
146  char expected = kUnset;
147  if(m_SCTBL.compare_exchange_strong(expected, kSetting)) {
149  m_SCTBL.store(kSet);
151  }
152  return tmp;
153 }
TrajectoryStateClosestToBeamLine trajectoryStateClosestToBeamLine
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
std::atomic< char > m_SCTBL
FreeTrajectoryState initialFTS
TrajectoryStateOnSurface TrackTransientTrack::stateOnSurface ( const GlobalPoint point) const
virtual

The TSOS at any point. The initial state will be used for the propagation.

Implements reco::BasicTransientTrack.

Definition at line 134 of file TrackTransientTrack.cc.

References TransverseImpactPointExtrapolator::extrapolate(), initialFTS, and theField.

135 {
137  return tipe.extrapolate(initialFTS, point);
138 }
const MagneticField * theField
FreeTrajectoryState initialFTS
const Track& reco::TrackTransientTrack::track ( ) const
inlinevirtual

Implements reco::BasicTransientTrack.

Definition at line 72 of file TrackTransientTrack.h.

72 {return *this;}
TrackBaseRef reco::TrackTransientTrack::trackBaseRef ( ) const
inlinevirtual

Implements reco::BasicTransientTrack.

Definition at line 66 of file TrackTransientTrack.h.

References tkr_.

Referenced by VertexFitterResult::fill().

66 {return TrackBaseRef(tkr_);}
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:32
TrajectoryStateClosestToPoint reco::TrackTransientTrack::trajectoryStateClosestToPoint ( const GlobalPoint point) const
inlinevirtual

Implements reco::BasicTransientTrack.

Definition at line 47 of file TrackTransientTrack.h.

References builder, and initialFTS.

48  {return builder(initialFTS, point);}
TSCPBuilderNoMaterial builder
FreeTrajectoryState initialFTS

Member Data Documentation

TSCPBuilderNoMaterial reco::TrackTransientTrack::builder
private

Definition at line 93 of file TrackTransientTrack.h.

Referenced by impactPointTSCP(), and trajectoryStateClosestToPoint().

FreeTrajectoryState reco::TrackTransientTrack::initialFTS
private
TrajectoryStateClosestToPoint reco::TrackTransientTrack::initialTSCP
mutableprivate

Definition at line 86 of file TrackTransientTrack.h.

Referenced by impactPointTSCP(), and TrackTransientTrack().

TrajectoryStateOnSurface reco::TrackTransientTrack::initialTSOS
mutableprivate

Definition at line 85 of file TrackTransientTrack.h.

Referenced by impactPointState(), and TrackTransientTrack().

std::atomic<char> reco::TrackTransientTrack::m_SCTBL
mutableprivate

Definition at line 91 of file TrackTransientTrack.h.

Referenced by setBeamSpot(), and stateAtBeamLine().

std::atomic<char> reco::TrackTransientTrack::m_TSCP
mutableprivate

Definition at line 90 of file TrackTransientTrack.h.

Referenced by impactPointTSCP(), and TrackTransientTrack().

std::atomic<char> reco::TrackTransientTrack::m_TSOS
mutableprivate
reco::BeamSpot reco::TrackTransientTrack::theBeamSpot
private

Definition at line 95 of file TrackTransientTrack.h.

Referenced by setBeamSpot(), and stateAtBeamLine().

const MagneticField* reco::TrackTransientTrack::theField
private
edm::ESHandle<GlobalTrackingGeometry> reco::TrackTransientTrack::theTrackingGeometry
private
TrackRef reco::TrackTransientTrack::tkr_
private

Definition at line 78 of file TrackTransientTrack.h.

Referenced by persistentTrackRef(), and trackBaseRef().

TrajectoryStateClosestToBeamLine reco::TrackTransientTrack::trajectoryStateClosestToBeamLine
mutableprivate

Definition at line 87 of file TrackTransientTrack.h.

Referenced by stateAtBeamLine().