FairRoot/PandaRoot
|
Abstract base class for a track representation. More...
#include <AbsTrackRep.h>
Public Member Functions | |
AbsTrackRep () | |
AbsTrackRep (int pdgCode, char propDir=0) | |
virtual | ~AbsTrackRep () |
virtual AbsTrackRep * | clone () const =0 |
Clone the trackRep. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
double | extrapolateToMeasurement (StateOnPlane &state, const AbsMeasurement *measurement, bool stopAtBoundary=false, bool calcJacobianNoise=false) const |
extrapolate to an AbsMeasurement More... | |
virtual unsigned int | getDim () const =0 |
Get the dimension of the state vector used by the track representation. More... | |
virtual TVector3 | getPos (const StateOnPlane &state) const =0 |
Get the cartesian position of a state. More... | |
virtual TVector3 | getMom (const StateOnPlane &state) const =0 |
Get the cartesian momentum vector of a state. More... | |
TVector3 | getDir (const StateOnPlane &state) const |
Get the direction vector of a state. More... | |
virtual void | getPosMom (const StateOnPlane &state, TVector3 &pos, TVector3 &mom) const =0 |
Get cartesian position and momentum vector of a state. More... | |
void | getPosDir (const StateOnPlane &state, TVector3 &pos, TVector3 &dir) const |
Get cartesian position and direction vector of a state. More... | |
virtual TVectorD | get6DState (const StateOnPlane &state) const |
Get the 6D state vector (x, y, z, p_x, p_y, p_z). More... | |
virtual TMatrixDSym | get6DCov (const MeasuredStateOnPlane &state) const =0 |
Get the 6D covariance. More... | |
virtual void | getPosMomCov (const MeasuredStateOnPlane &state, TVector3 &pos, TVector3 &mom, TMatrixDSym &cov) const =0 |
Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance. More... | |
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. More... | |
virtual double | getMomMag (const StateOnPlane &state) const =0 |
get the magnitude of the momentum in GeV. More... | |
virtual double | getMomVar (const MeasuredStateOnPlane &state) const =0 |
get the variance of the absolute value of the momentum . More... | |
int | getPDG () const |
Get the pdg code. More... | |
double | getPDGCharge () const |
Get the charge of the particle of the pdg code. More... | |
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). More... | |
virtual double | getQop (const StateOnPlane &state) const =0 |
Get charge over momentum. More... | |
double | getMass (const StateOnPlane &state) const |
Get tha particle mass in GeV/c^2. More... | |
char | getPropDir () const |
Get propagation direction. (-1, 0, 1) -> (backward, auto, forward). More... | |
virtual void | getForwardJacobianAndNoise (TMatrixD &jacobian, TMatrixDSym &noise, TVectorD &deltaState) const =0 |
Get the jacobian and noise matrix of the last extrapolation. More... | |
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. More... | |
virtual std::vector < genfit::MatStep > | getSteps () const =0 |
Get stepsizes and material properties of crossed materials of the last extrapolation. More... | |
virtual double | getRadiationLenght () const =0 |
Get the accumulated X/X0 (path / radiation length) of the material crossed in the last extrapolation. More... | |
virtual double | getTime (const StateOnPlane &) const =0 |
Get the time corresponding to the StateOnPlane. Extrapolation. More... | |
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. More... | |
bool | switchPDGSign () |
try to multiply pdg code with -1. (Switch from particle to anti-particle and vice versa). More... | |
virtual void | setPosMom (StateOnPlane &state, const TVector3 &pos, const TVector3 &mom) const =0 |
Set position and momentum of state. More... | |
virtual void | setPosMom (StateOnPlane &state, const TVectorD &state6) const =0 |
Set position and momentum of state. More... | |
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. More... | |
virtual void | setPosMomCov (MeasuredStateOnPlane &state, const TVector3 &pos, const TVector3 &mom, const TMatrixDSym &cov6x6) const =0 |
Set position, momentum and covariance of state. More... | |
virtual void | setPosMomCov (MeasuredStateOnPlane &state, const TVectorD &state6, const TMatrixDSym &cov6x6) const =0 |
Set position, momentum and covariance of state. More... | |
virtual void | setChargeSign (StateOnPlane &state, double charge) const =0 |
Set the sign of the charge according to charge. More... | |
virtual void | setQop (StateOnPlane &state, double qop) const =0 |
Set charge/momentum. More... | |
virtual void | setTime (StateOnPlane &state, double time) const =0 |
Set time at which the state was defined. More... | |
void | setPropDir (int dir) |
Set propagation direction. (-1, 0, 1) -> (backward, auto, forward). More... | |
void | switchPropDir () |
Switch propagation direction. Has no effect if propDir_ is set to 0. More... | |
virtual bool | isSameType (const AbsTrackRep *other)=0 |
check if other is of same type (e.g. RKTrackRep). More... | |
virtual bool | isSame (const AbsTrackRep *other)=0 |
check if other is of same type (e.g. RKTrackRep) and has same pdg code. More... | |
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! More... | |
AbsTrackRep & | operator= (const AbsTrackRep &) |
protect from calling assignment operator from outside the class. Use clone() instead! More... | |
Protected Attributes | |
int | pdgCode_ |
Particle code. More... | |
char | propDir_ |
propagation direction (-1, 0, 1) -> (backward, auto, forward) More... | |
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.
genfit::AbsTrackRep::AbsTrackRep | ( | ) |
genfit::AbsTrackRep::AbsTrackRep | ( | int | pdgCode, |
char | propDir = 0 |
||
) |
|
inlinevirtual |
Definition 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!
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.
|
pure virtual |
Clone the trackRep.
Implemented in genfit::RKTrackRep.
|
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 genfit::StateOnPlane::extrapolateBy().
|
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 genfit::StateOnPlane::extrapolateToCylinder().
|
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 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.
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
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 genfit::StateOnPlane::extrapolateToPlane().
|
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 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 genfit::StateOnPlane::extrapolateToSphere().
|
pure virtual |
Get the 6D covariance.
Implemented in genfit::RKTrackRep.
Referenced by genfit::MeasuredStateOnPlane::get6DCov().
|
virtual |
Get the 6D state vector (x, y, z, p_x, p_y, p_z).
Referenced by genfit::StateOnPlane::get6DState().
|
virtual |
Translates MeasuredStateOnPlane into 6D state vector (x, y, z, p_x, p_y, p_z) and 6x6 covariance.
Referenced by genfit::MeasuredStateOnPlane::get6DStateCov().
|
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.
|
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.
Referenced by genfit::StateOnPlane::getCharge().
|
pure virtual |
Get the dimension of the state vector used by the track representation.
Implemented in genfit::RKTrackRep.
Referenced by genfit::MeasuredStateOnPlane::MeasuredStateOnPlane(), and genfit::StateOnPlane::StateOnPlane().
|
inline |
Get the direction vector of a state.
Definition at line 230 of file AbsTrackRep.h.
References getMom().
Referenced by genfit::StateOnPlane::getDir().
|
pure virtual |
Get the jacobian and noise matrix of the last extrapolation.
Implemented in genfit::RKTrackRep.
double genfit::AbsTrackRep::getMass | ( | const StateOnPlane & | state | ) | const |
Get tha particle mass in GeV/c^2.
Referenced by genfit::StateOnPlane::getMass().
|
pure virtual |
Get the cartesian momentum vector of a state.
Implemented in genfit::RKTrackRep.
Referenced by extrapolateToLine(), getDir(), and genfit::StateOnPlane::getMom().
|
pure virtual |
get the magnitude of the momentum in GeV.
Implemented in genfit::RKTrackRep.
Referenced by genfit::StateOnPlane::getMomMag().
|
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 256 of file AbsTrackRep.h.
References pdgCode_.
Referenced by genfit::StateOnPlane::getPDG().
double genfit::AbsTrackRep::getPDGCharge | ( | ) | const |
Get the charge of the particle of the pdg code.
|
pure virtual |
Get the cartesian position of a state.
Implemented in genfit::RKTrackRep.
Referenced by extrapolateToLine(), and genfit::StateOnPlane::getPos().
|
inline |
Get cartesian position and direction vector of a state.
Definition at line 236 of file AbsTrackRep.h.
References getPosMom().
Referenced by genfit::StateOnPlane::getPosDir().
|
pure virtual |
Get cartesian position and momentum vector of a state.
Implemented in genfit::RKTrackRep.
Referenced by getPosDir(), and genfit::StateOnPlane::getPosMom().
|
pure virtual |
Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance.
Implemented in genfit::RKTrackRep.
Referenced by genfit::MeasuredStateOnPlane::getPosMomCov().
|
inline |
Get propagation direction. (-1, 0, 1) -> (backward, auto, forward).
Definition at line 272 of file AbsTrackRep.h.
References propDir_.
|
pure virtual |
Get charge over momentum.
Implemented in genfit::RKTrackRep.
Referenced by genfit::StateOnPlane::getQop().
|
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.
|
pure virtual |
Get the time corresponding to the StateOnPlane. Extrapolation.
Implemented in genfit::RKTrackRep.
Referenced by genfit::StateOnPlane::getTime().
|
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 |
|
pure virtual |
Set the sign of the charge according to charge.
Implemented in genfit::RKTrackRep.
Referenced by genfit::StateOnPlane::setChargeSign().
|
inlinevirtual |
|
pure virtual |
Set position and momentum of state.
Implemented in genfit::RKTrackRep.
Referenced by genfit::StateOnPlane::setPosMom().
|
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 genfit::MeasuredStateOnPlane::setPosMomCov().
|
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 320 of file AbsTrackRep.h.
References propDir_.
|
pure virtual |
Set charge/momentum.
Implemented in genfit::RKTrackRep.
Referenced by genfit::StateOnPlane::setQop().
|
pure virtual |
Set time at which the state was defined.
Implemented in genfit::RKTrackRep.
Referenced by genfit::StateOnPlane::setTime().
bool genfit::AbsTrackRep::switchPDGSign | ( | ) |
try to multiply pdg code with -1. (Switch from particle to anti-particle and vice versa).
|
inline |
Switch propagation direction. Has no effect if propDir_ is set to 0.
Definition at line 327 of file AbsTrackRep.h.
References propDir_.
|
protected |
Definition at line 352 of file AbsTrackRep.h.
Referenced by setDebugLvl().
|
protected |
|
protected |
propagation direction (-1, 0, 1) -> (backward, auto, forward)
Definition at line 350 of file AbsTrackRep.h.
Referenced by getPropDir(), setPropDir(), and switchPropDir().