Analysis Software
Documentation for sPHENIX simulation software
|
Abstract base class for a track representation. More...
#include <GenFit/blob/master/core/include/AbsTrackRep.h>
Public Member Functions | |
AbsTrackRep () | |
AbsTrackRep (int pdgCode, char propDir=0) | |
virtual | ~AbsTrackRep () |
virtual AbsTrackRep * | clone () const =0 |
Clone the trackRep. | |
virtual double | extrapolateToPlane (StateOnPlane &state, const genfit::SharedPlanePtr &plane, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0 |
Extrapolates the state to plane, and returns the extrapolation length and, via reference, the extrapolated state. | |
virtual double | extrapolateToLine (StateOnPlane &state, const TVector3 &linePoint, const TVector3 &lineDirection, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0 |
Extrapolates the state to the POCA to a line, and returns the extrapolation length and, via reference, the extrapolated state. | |
virtual double | extrapolateToLine (StateOnPlane &state, const TVector3 &point1, const TVector3 &point2, TVector3 &poca, TVector3 &dirInPoca, TVector3 &poca_onwire, bool stopAtBoundary=false, bool calcJacobianNoise=false) const |
Resembles the interface of GFAbsTrackRep in old versions of genfit. | |
virtual double | extrapolateToPoint (StateOnPlane &state, const TVector3 &point, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0 |
Extrapolates the state to the POCA to a point, and returns the extrapolation length and, via reference, the extrapolated state. | |
virtual double | extrapolateToPoint (StateOnPlane &state, const TVector3 &point, const TMatrixDSym &G, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0 |
Extrapolates the state to the POCA to a point in the metric of G, and returns the extrapolation length and, via reference, the extrapolated state. | |
virtual double | extrapolateToCylinder (StateOnPlane &state, double radius, const TVector3 &linePoint=TVector3(0., 0., 0.), const TVector3 &lineDirection=TVector3(0., 0., 1.), bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0 |
Extrapolates the state to the cylinder surface, and returns the extrapolation length and, via reference, the extrapolated state. | |
virtual double | extrapolateToCone (StateOnPlane &state, double radius, const TVector3 &linePoint=TVector3(0., 0., 0.), const TVector3 &lineDirection=TVector3(0., 0., 1.), bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0 |
Extrapolates the state to the cone surface, and returns the extrapolation length and, via reference, the extrapolated state. | |
virtual double | extrapolateToSphere (StateOnPlane &state, double radius, const TVector3 &point=TVector3(0., 0., 0.), bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0 |
Extrapolates the state to the sphere surface, and returns the extrapolation length and, via reference, the extrapolated state. | |
virtual double | extrapolateBy (StateOnPlane &state, double step, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0 |
Extrapolates the state by step (cm) and returns the extrapolation length and, via reference, the extrapolated state. | |
double | extrapolateToMeasurement (StateOnPlane &state, const AbsMeasurement *measurement, bool stopAtBoundary=false, bool calcJacobianNoise=false) const |
extrapolate to an AbsMeasurement | |
virtual unsigned int | getDim () const =0 |
Get the dimension of the state vector used by the track representation. | |
virtual TVector3 | getPos (const StateOnPlane &state) const =0 |
Get the cartesian position of a state. | |
virtual TVector3 | getMom (const StateOnPlane &state) const =0 |
Get the cartesian momentum vector of a state. | |
TVector3 | getDir (const StateOnPlane &state) const |
Get the direction vector of a state. | |
virtual void | getPosMom (const StateOnPlane &state, TVector3 &pos, TVector3 &mom) const =0 |
Get cartesian position and momentum vector of a state. | |
void | getPosDir (const StateOnPlane &state, TVector3 &pos, TVector3 &dir) const |
Get cartesian position and direction vector of a state. | |
virtual TVectorD | get6DState (const StateOnPlane &state) const |
Get the 6D state vector (x, y, z, p_x, p_y, p_z). | |
virtual TMatrixDSym | get6DCov (const MeasuredStateOnPlane &state) const =0 |
Get the 6D covariance. | |
virtual void | getPosMomCov (const MeasuredStateOnPlane &state, TVector3 &pos, TVector3 &mom, TMatrixDSym &cov) const =0 |
Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance. | |
virtual void | get6DStateCov (const MeasuredStateOnPlane &state, TVectorD &stateVec, TMatrixDSym &cov) const |
Translates MeasuredStateOnPlane into 6D state vector (x, y, z, p_x, p_y, p_z) and 6x6 covariance. | |
virtual double | getMomMag (const StateOnPlane &state) const =0 |
get the magnitude of the momentum in GeV. | |
virtual double | getMomVar (const MeasuredStateOnPlane &state) const =0 |
get the variance of the absolute value of the momentum . | |
int | getPDG () const |
Get the pdg code. | |
double | getPDGCharge () const |
Get the charge of the particle of the pdg code. | |
virtual double | getCharge (const StateOnPlane &state) const =0 |
Get the (fitted) charge of a state. This is not always equal the pdg charge (e.g. if the charge sign was flipped during the fit). | |
virtual double | getQop (const StateOnPlane &state) const =0 |
Get charge over momentum. | |
double | getMass (const StateOnPlane &state) const |
Get tha particle mass in GeV/c^2. | |
char | getPropDir () const |
Get propagation direction. (-1, 0, 1) -> (backward, auto, forward). | |
virtual void | getForwardJacobianAndNoise (TMatrixD &jacobian, TMatrixDSym &noise, TVectorD &deltaState) const =0 |
Get the jacobian and noise matrix of the last extrapolation. | |
virtual void | getBackwardJacobianAndNoise (TMatrixD &jacobian, TMatrixDSym &noise, TVectorD &deltaState) const =0 |
Get the jacobian and noise matrix of the last extrapolation if it would have been done in opposite direction. | |
virtual std::vector < genfit::MatStep > | getSteps () const =0 |
Get stepsizes and material properties of crossed materials of the last extrapolation. | |
virtual double | getRadiationLenght () const =0 |
Get the accumulated X/X0 (path / radiation length) of the material crossed in the last extrapolation. | |
virtual double | getTime (const StateOnPlane &) const =0 |
Get the time corresponding to the StateOnPlane. Extrapolation. | |
void | calcJacobianNumerically (const genfit::StateOnPlane &origState, const genfit::SharedPlanePtr destPlane, TMatrixD &jacobian) const |
Calculate Jacobian of transportation numerically. Slow but accurate. Can be used to validate (semi)analytic calculations. | |
bool | switchPDGSign () |
try to multiply pdg code with -1. (Switch from particle to anti-particle and vice versa). | |
virtual void | setPosMom (StateOnPlane &state, const TVector3 &pos, const TVector3 &mom) const =0 |
Set position and momentum of state. | |
virtual void | setPosMom (StateOnPlane &state, const TVectorD &state6) const =0 |
Set position and momentum of state. | |
virtual void | setPosMomErr (MeasuredStateOnPlane &state, const TVector3 &pos, const TVector3 &mom, const TVector3 &posErr, const TVector3 &momErr) const =0 |
Set position and momentum and error of state. | |
virtual void | setPosMomCov (MeasuredStateOnPlane &state, const TVector3 &pos, const TVector3 &mom, const TMatrixDSym &cov6x6) const =0 |
Set position, momentum and covariance of state. | |
virtual void | setPosMomCov (MeasuredStateOnPlane &state, const TVectorD &state6, const TMatrixDSym &cov6x6) const =0 |
Set position, momentum and covariance of state. | |
virtual void | setChargeSign (StateOnPlane &state, double charge) const =0 |
Set the sign of the charge according to charge. | |
virtual void | setQop (StateOnPlane &state, double qop) const =0 |
Set charge/momentum. | |
virtual void | setTime (StateOnPlane &state, double time) const =0 |
Set time at which the state was defined. | |
void | setPropDir (int dir) |
Set propagation direction. (-1, 0, 1) -> (backward, auto, forward). | |
void | switchPropDir () |
Switch propagation direction. Has no effect if propDir_ is set to 0. | |
virtual bool | isSameType (const AbsTrackRep *other)=0 |
check if other is of same type (e.g. RKTrackRep). | |
virtual bool | isSame (const AbsTrackRep *other)=0 |
check if other is of same type (e.g. RKTrackRep) and has same pdg code. | |
virtual void | setDebugLvl (unsigned int lvl=1) |
virtual void | Print (const Option_t *="") const |
Protected Member Functions | |
AbsTrackRep (const AbsTrackRep &) | |
protect from calling copy c'tor from outside the class. Use clone() if you want a copy! | |
AbsTrackRep & | operator= (const AbsTrackRep &) |
protect from calling assignment operator from outside the class. Use clone() instead! | |
Protected Attributes | |
int | pdgCode_ |
Particle code. | |
char | propDir_ |
propagation direction (-1, 0, 1) -> (backward, auto, forward) | |
unsigned int | debugLvl_ |
Abstract base class for a track representation.
Provides functionality to extrapolate a StateOnPlane to another DetPlane, to the POCA to a line or a point, or a cylinder or sphere. Defines a set of parameters describing the track. StateOnPlane objects are always defined with a track parameterization of a specific AbsTrackRep. The AbsTrackRep provides functionality to translate from the internal representation of a state into cartesian position and momentum (and covariance) and vice versa.
Definition at line 66 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 66 of file AbsTrackRep.h
genfit::AbsTrackRep::AbsTrackRep | ( | ) |
Definition at line 31 of file AbsTrackRep.cc.
View newest version in sPHENIX GitHub at line 31 of file AbsTrackRep.cc
genfit::AbsTrackRep::AbsTrackRep | ( | int | pdgCode, |
char | propDir = 0 |
||
) |
Definition at line 37 of file AbsTrackRep.cc.
View newest version in sPHENIX GitHub at line 37 of file AbsTrackRep.cc
|
inlinevirtual |
Definition at line 73 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 73 of file AbsTrackRep.h
|
protected |
protect from calling copy c'tor from outside the class. Use clone() if you want a copy!
Definition at line 43 of file AbsTrackRep.cc.
View newest version in sPHENIX GitHub at line 43 of file AbsTrackRep.cc
void genfit::AbsTrackRep::calcJacobianNumerically | ( | const genfit::StateOnPlane & | origState, |
const genfit::SharedPlanePtr | destPlane, | ||
TMatrixD & | jacobian | ||
) | const |
Calculate Jacobian of transportation numerically. Slow but accurate. Can be used to validate (semi)analytic calculations.
Definition at line 105 of file AbsTrackRep.cc.
View newest version in sPHENIX GitHub at line 105 of file AbsTrackRep.cc
References extrapolateToPlane(), getDim(), genfit::StateOnPlane::getState(), i, and j.
Referenced by compareForthBackJacNoise().
|
pure virtual |
Clone the trackRep.
Implemented in genfit::RKTrackRep.
Referenced by main(), and genfit::Track::Track().
|
pure virtual |
Extrapolates the state by step (cm) and returns the extrapolation length and, via reference, the extrapolated state.
If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.
If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.
Implemented in genfit::RKTrackRep.
Referenced by checkExtrapolateBy(), genfit::GblFitter::constructGblInfo(), genfit::StateOnPlane::extrapolateBy(), and genfit::GFGbl::processTrackWithRep().
|
pure virtual |
Extrapolates the state to the cone surface, and returns the extrapolation length and, via reference, the extrapolated state.
If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.
If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.
Implemented in genfit::RKTrackRep.
Referenced by genfit::StateOnPlane::extrapolateToCone().
|
pure virtual |
Extrapolates the state to the cylinder surface, and returns the extrapolation length and, via reference, the extrapolated state.
If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.
If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.
Implemented in genfit::RKTrackRep.
Referenced by checkExtrapolateToCylinder(), PHGenFit::extrapolateToCylinder(), PHGenFit::Track::extrapolateToCylinder(), genfit::StateOnPlane::extrapolateToCylinder(), and GenFitTrackProp::fill_tree().
|
pure virtual |
Extrapolates the state to the POCA to a line, and returns the extrapolation length and, via reference, the extrapolated state.
If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.
If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.
Implemented in genfit::RKTrackRep.
Referenced by checkExtrapolateToLine(), genfit::ProlateSpacepointMeasurement::constructPlane(), genfit::WirePointMeasurement::constructPlane(), genfit::WireMeasurement::constructPlane(), genfit::WireMeasurementNew::constructPlane(), PHGenFit::Track::extrapolateToLine(), genfit::StateOnPlane::extrapolateToLine(), and extrapolateToLine().
|
inlinevirtual |
Resembles the interface of GFAbsTrackRep in old versions of genfit.
This interface to extrapolateToLine is intended to resemble the interface of GFAbsTrackRep in old versions of genfit and is implemented by default via the preceding function.
If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.
If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.
Definition at line 120 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 120 of file AbsTrackRep.h
References extrapolateToLine(), getMom(), and getPos().
double genfit::AbsTrackRep::extrapolateToMeasurement | ( | StateOnPlane & | state, |
const AbsMeasurement * | measurement, | ||
bool | stopAtBoundary = false , |
||
bool | calcJacobianNoise = false |
||
) | const |
extrapolate to an AbsMeasurement
Definition at line 50 of file AbsTrackRep.cc.
View newest version in sPHENIX GitHub at line 50 of file AbsTrackRep.cc
References genfit::AbsMeasurement::constructPlane(), and extrapolateToPlane().
Referenced by genfit::StateOnPlane::extrapolateToMeasurement().
|
pure virtual |
Extrapolates the state to plane, and returns the extrapolation length and, via reference, the extrapolated state.
If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.
If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.
Implemented in genfit::RKTrackRep.
Referenced by calcJacobianNumerically(), checkErrorPropagation(), checkStopAtBoundary(), compareForthBackExtrapolation(), compareForthBackJacNoise(), genfit::GblFitter::constructGblInfo(), extrapolateToMeasurement(), PHGenFit::Track::extrapolateToPlane(), genfit::StateOnPlane::extrapolateToPlane(), GenFitTrackProp::fill_tree(), genfit::Track::getTrackLen(), main(), genfit::EventDisplay::makeLines(), genfit::KalmanFitterRefTrack::prepareTrack(), genfit::KalmanFitter::processTrackPoint(), genfit::GFGbl::processTrackWithRep(), genfit::GblFitterInfo::recalculateJacobian(), and PHGenFit::Track::updateOneMeasurementKalman().
|
pure virtual |
Extrapolates the state to the POCA to a point, and returns the extrapolation length and, via reference, the extrapolated state.
If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.
If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.
Implemented in genfit::RKTrackRep.
Referenced by checkExtrapolateToPoint(), PHGenFit::Track::extrapolateToPoint(), and genfit::StateOnPlane::extrapolateToPoint().
|
pure virtual |
Extrapolates the state to the POCA to a point in the metric of G, and returns the extrapolation length and, via reference, the extrapolated state.
If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.
If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.
Implemented in genfit::RKTrackRep.
|
pure virtual |
Extrapolates the state to the sphere surface, and returns the extrapolation length and, via reference, the extrapolated state.
If stopAtBoundary is true, the extrapolation stops as soon as a material boundary is encountered.
If state has a covariance, jacobian and noise matrices will be calculated and the covariance will be propagated. If state has no covariance, jacobian and noise will only be calculated if calcJacobianNoise == true.
Implemented in genfit::RKTrackRep.
Referenced by checkExtrapolateToSphere(), and genfit::StateOnPlane::extrapolateToSphere().
|
pure virtual |
Get the 6D covariance.
Implemented in genfit::RKTrackRep.
Referenced by genfit::MeasuredStateOnPlane::get6DCov(), TrackProjectorPid::project_track(), and TrackProjectorPlaneECAL::project_track().
|
virtual |
Get the 6D state vector (x, y, z, p_x, p_y, p_z).
Definition at line 59 of file AbsTrackRep.cc.
View newest version in sPHENIX GitHub at line 59 of file AbsTrackRep.cc
References getPosMom(), and Acts::Test::pos.
Referenced by genfit::StateOnPlane::get6DState().
|
virtual |
Translates MeasuredStateOnPlane into 6D state vector (x, y, z, p_x, p_y, p_z) and 6x6 covariance.
Definition at line 77 of file AbsTrackRep.cc.
View newest version in sPHENIX GitHub at line 77 of file AbsTrackRep.cc
References getPosMomCov(), and Acts::Test::pos.
Referenced by genfit::MeasuredStateOnPlane::get6DStateCov(), and main().
|
pure virtual |
Get the jacobian and noise matrix of the last extrapolation if it would have been done in opposite direction.
Implemented in genfit::RKTrackRep.
Referenced by compareForthBackJacNoise(), and genfit::KalmanFitterRefTrack::prepareTrack().
|
pure virtual |
Get the (fitted) charge of a state. This is not always equal the pdg charge (e.g. if the charge sign was flipped during the fit).
Implemented in genfit::RKTrackRep, and genfit::MplTrackRep.
Referenced by genfit::GblFitter::constructGblInfo(), PHGenFit::Track::get_charge(), genfit::StateOnPlane::getCharge(), genfit::KalmanFitterRefTrack::prepareTrack(), genfit::KalmanFitterRefTrack::processTrackPoint(), genfit::KalmanFitterRefTrack::processTrackPointSqrt(), genfit::KalmanFitterRefTrack::processTrackWithRep(), and genfit::GFGbl::processTrackWithRep().
|
pure virtual |
Get the dimension of the state vector used by the track representation.
Implemented in genfit::RKTrackRep.
Referenced by calcJacobianNumerically(), genfit::KalmanFitterInfo::checkConsistency(), genfit::GblFitter::constructGblInfo(), genfit::EventDisplay::drawEvent(), genfit::KalmanFitterRefTrack::fitTrack(), genfit::KalmanFitter::fitTrack(), genfit::FullMeasurement::FullMeasurement(), genfit::GblFitterInfo::GblFitterInfo(), genfit::AbsKalmanFitter::getChiSquNdf(), genfit::MeasuredStateOnPlane::MeasuredStateOnPlane(), genfit::KalmanFitterRefTrack::prepareTrack(), genfit::KalmanFitterRefTrack::processTrackPoint(), genfit::KalmanFitterRefTrack::processTrackPointSqrt(), genfit::GFGbl::processTrackWithRep(), and genfit::StateOnPlane::StateOnPlane().
|
inline |
Get the direction vector of a state.
Definition at line 246 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 246 of file AbsTrackRep.h
References getMom().
Referenced by genfit::GblFitter::constructGblInfo(), genfit::StateOnPlane::getDir(), and genfit::GFGbl::processTrackWithRep().
|
pure virtual |
Get the jacobian and noise matrix of the last extrapolation.
Implemented in genfit::RKTrackRep.
Referenced by compareForthBackJacNoise(), genfit::GblFitter::constructGblInfo(), genfit::GblFitterInfo::GblFitterInfo(), genfit::KalmanFitterRefTrack::prepareTrack(), genfit::GFGbl::processTrackWithRep(), and genfit::GblFitterInfo::recalculateJacobian().
double genfit::AbsTrackRep::getMass | ( | const StateOnPlane & | state | ) | const |
Get tha particle mass in GeV/c^2.
Definition at line 100 of file AbsTrackRep.cc.
View newest version in sPHENIX GitHub at line 100 of file AbsTrackRep.cc
References pdgCode_.
Referenced by genfit::GblFitter::constructGblInfo(), genfit::RKTrackRep::extrapolateBy(), genfit::RKTrackRep::extrapolateToCone(), genfit::RKTrackRep::extrapolateToCylinder(), genfit::RKTrackRep::extrapolateToLine(), genfit::RKTrackRep::extrapolateToPlane(), genfit::RKTrackRep::extrapolateToSphere(), genfit::RKTrackRep::extrapToPoint(), genfit::StateOnPlane::getMass(), and genfit::GFGbl::processTrackWithRep().
|
pure virtual |
Get the cartesian momentum vector of a state.
Implemented in genfit::RKTrackRep.
Referenced by checkExtrapolateToLine(), checkExtrapolateToPoint(), checkSetGetPosMom(), compareForthBackExtrapolation(), genfit::ProlateSpacepointMeasurement::constructPlane(), genfit::WirePointMeasurement::constructPlane(), genfit::WireMeasurement::constructPlane(), genfit::WireMeasurementNew::constructPlane(), extrapolateToLine(), getDir(), genfit::StateOnPlane::getMom(), PHGenFitTrackProjection::process_event(), TrackProjectorPid::project_track(), and TrackProjectorPlaneECAL::project_track().
|
pure virtual |
get the magnitude of the momentum in GeV.
Implemented in genfit::RKTrackRep.
Referenced by genfit::GblFitter::constructGblInfo(), genfit::StateOnPlane::getMomMag(), and genfit::GFGbl::processTrackWithRep().
|
pure virtual |
get the variance of the absolute value of the momentum .
Implemented in genfit::RKTrackRep.
Referenced by genfit::MeasuredStateOnPlane::getMomVar().
|
inline |
Get the pdg code.
Definition at line 272 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 272 of file AbsTrackRep.h
References pdgCode_.
Referenced by genfit::GFRaveTrackParameters::getPdg(), genfit::StateOnPlane::getPDG(), genfit::RKTrackRep::isSame(), and genfit::GFGbl::processTrackWithRep().
double genfit::AbsTrackRep::getPDGCharge | ( | ) | const |
Get the charge of the particle of the pdg code.
Definition at line 93 of file AbsTrackRep.cc.
View newest version in sPHENIX GitHub at line 93 of file AbsTrackRep.cc
References assert, particle, and pdgCode_.
Referenced by genfit::MplTrackRep::getCharge(), and genfit::RKTrackRep::getCharge().
|
pure virtual |
Get the cartesian position of a state.
Implemented in genfit::RKTrackRep.
Referenced by checkExtrapolateToCylinder(), checkExtrapolateToSphere(), checkSetGetPosMom(), checkStopAtBoundary(), compareForthBackExtrapolation(), genfit::WireMeasurement::constructPlane(), genfit::WireMeasurementNew::constructPlane(), extrapolateToLine(), genfit::StateOnPlane::getPos(), PHGenFitTrackProjection::process_event(), TrackProjectorPid::project_track(), and TrackProjectorPlaneECAL::project_track().
|
inline |
Get cartesian position and direction vector of a state.
Definition at line 252 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 252 of file AbsTrackRep.h
References getPosMom().
Referenced by genfit::StateOnPlane::getPosDir(), and genfit::EventDisplay::makeLines().
|
pure virtual |
Get cartesian position and momentum vector of a state.
Implemented in genfit::RKTrackRep.
Referenced by get6DState(), getPosDir(), genfit::StateOnPlane::getPosMom(), genfit::KalmanFitterRefTrack::prepareTrack(), genfit::StateOnPlane::Print(), genfit::KalmanFitterRefTrack::processTrackPoint(), and genfit::KalmanFitterRefTrack::processTrackPointSqrt().
|
pure virtual |
Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance.
Implemented in genfit::RKTrackRep.
Referenced by get6DStateCov(), genfit::MeasuredStateOnPlane::getPosMomCov(), genfit::EventDisplay::makeLines(), genfit::MeasuredStateOnPlane::Print(), genfit::KalmanFitter::processTrackPartially(), and genfit::KalmanFitter::processTrackWithRep().
|
inline |
Get propagation direction. (-1, 0, 1) -> (backward, auto, forward).
Definition at line 288 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 288 of file AbsTrackRep.h
References propDir_.
|
pure virtual |
Get charge over momentum.
Implemented in genfit::RKTrackRep.
Referenced by genfit::StateOnPlane::getQop(), genfit::KalmanFitterRefTrack::prepareTrack(), genfit::KalmanFitter::processTrackPartially(), and genfit::KalmanFitter::processTrackWithRep().
|
pure virtual |
Get the accumulated X/X0 (path / radiation length) of the material crossed in the last extrapolation.
Implemented in genfit::RKTrackRep.
|
pure virtual |
Get stepsizes and material properties of crossed materials of the last extrapolation.
Implemented in genfit::RKTrackRep.
Referenced by genfit::GblFitter::constructGblInfo(), and genfit::GFGbl::processTrackWithRep().
|
pure virtual |
Get the time corresponding to the StateOnPlane. Extrapolation.
Implemented in genfit::RKTrackRep.
Referenced by genfit::StateOnPlane::getTime(), and genfit::Track::getTOF().
|
pure virtual |
check if other is of same type (e.g. RKTrackRep) and has same pdg code.
Implemented in genfit::RKTrackRep.
|
pure virtual |
check if other is of same type (e.g. RKTrackRep).
Implemented in genfit::RKTrackRep.
|
protected |
protect from calling assignment operator from outside the class. Use clone() instead!
|
virtual |
Definition at line 194 of file AbsTrackRep.cc.
View newest version in sPHENIX GitHub at line 194 of file AbsTrackRep.cc
References pdgCode_, genfit::printOut, and propDir_.
|
pure virtual |
Set the sign of the charge according to charge.
Implemented in genfit::RKTrackRep.
Referenced by genfit::StateOnPlane::setChargeSign().
|
inlinevirtual |
Definition at line 351 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 351 of file AbsTrackRep.h
References debugLvl_.
|
pure virtual |
Set position and momentum of state.
Implemented in genfit::RKTrackRep.
Referenced by checkErrorPropagation(), checkExtrapolateBy(), checkExtrapolateToCylinder(), checkExtrapolateToLine(), checkExtrapolateToPoint(), checkExtrapolateToSphere(), checkSetGetPosMom(), checkStopAtBoundary(), compareForthBackExtrapolation(), compareForthBackJacNoise(), genfit::GblFitter::constructGblInfo(), genfit::KalmanFitterRefTrack::prepareTrack(), genfit::StateOnPlane::setPosMom(), and genfit::GblFitter::sortHits().
|
pure virtual |
Set position and momentum of state.
Implemented in genfit::RKTrackRep.
|
pure virtual |
Set position, momentum and covariance of state.
Implemented in genfit::RKTrackRep.
Referenced by PHGenFit::Track::extrapolateToCylinder(), main(), PHGenFitTrackProjection::process_event(), genfit::KalmanFitter::processTrackPartially(), genfit::KalmanFitterRefTrack::processTrackPoint(), genfit::KalmanFitterRefTrack::processTrackPointSqrt(), genfit::KalmanFitter::processTrackWithRep(), TrackProjectorPid::project_track(), TrackProjectorPlaneECAL::project_track(), genfit::MeasuredStateOnPlane::setPosMomCov(), and PHGenFit::Track::updateOneMeasurementKalman().
|
pure virtual |
Set position, momentum and covariance of state.
Implemented in genfit::RKTrackRep.
|
pure virtual |
Set position and momentum and error of state.
Implemented in genfit::RKTrackRep.
Referenced by genfit::MeasuredStateOnPlane::setPosMomErr().
|
inline |
Set propagation direction. (-1, 0, 1) -> (backward, auto, forward).
Definition at line 336 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 336 of file AbsTrackRep.h
References propDir_.
Referenced by genfit::EventDisplay::drawEvent().
|
pure virtual |
Set charge/momentum.
Implemented in genfit::RKTrackRep.
Referenced by genfit::KalmanFitterRefTrack::prepareTrack(), genfit::KalmanFitter::processTrackPartially(), genfit::KalmanFitter::processTrackWithRep(), and genfit::StateOnPlane::setQop().
|
pure virtual |
Set time at which the state was defined.
Implemented in genfit::RKTrackRep.
Referenced by genfit::GblFitter::constructGblInfo(), genfit::KalmanFitterRefTrack::prepareTrack(), genfit::KalmanFitter::processTrackPartially(), genfit::KalmanFitter::processTrackWithRep(), genfit::StateOnPlane::setTime(), and genfit::GblFitter::sortHits().
bool genfit::AbsTrackRep::switchPDGSign | ( | ) |
try to multiply pdg code with -1. (Switch from particle to anti-particle and vice versa).
Definition at line 183 of file AbsTrackRep.cc.
View newest version in sPHENIX GitHub at line 183 of file AbsTrackRep.cc
References particle, and pdgCode_.
Referenced by genfit::Track::switchPDGSigns().
|
inline |
Switch propagation direction. Has no effect if propDir_ is set to 0.
Definition at line 343 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 343 of file AbsTrackRep.h
References propDir_.
|
protected |
Definition at line 368 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 368 of file AbsTrackRep.h
Referenced by genfit::RKTrackRep::calcForwardJacobianAndNoise(), genfit::RKTrackRep::checkCache(), genfit::RKTrackRep::estimateStep(), genfit::RKTrackRep::Extrap(), genfit::RKTrackRep::extrapolateBy(), genfit::RKTrackRep::extrapolateToCone(), genfit::RKTrackRep::extrapolateToCylinder(), genfit::RKTrackRep::extrapolateToLine(), genfit::RKTrackRep::extrapolateToPlane(), genfit::RKTrackRep::extrapolateToSphere(), genfit::RKTrackRep::extrapToPoint(), genfit::RKTrackRep::getBackwardJacobianAndNoise(), genfit::RKTrackRep::getForwardJacobianAndNoise(), genfit::MplTrackRep::RKPropagate(), genfit::RKTrackRep::RKPropagate(), genfit::RKTrackRep::RKutta(), and setDebugLvl().
|
protected |
Particle code.
Definition at line 364 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 364 of file AbsTrackRep.h
Referenced by genfit::RKTrackRep::estimateStep(), genfit::RKTrackRep::Extrap(), getMass(), getPDG(), getPDGCharge(), Print(), and switchPDGSign().
|
protected |
propagation direction (-1, 0, 1) -> (backward, auto, forward)
Definition at line 366 of file AbsTrackRep.h.
View newest version in sPHENIX GitHub at line 366 of file AbsTrackRep.h
Referenced by genfit::RKTrackRep::estimateStep(), getPropDir(), Print(), setPropDir(), and switchPropDir().