SND@LHC Software
Loading...
Searching...
No Matches
genfit::AbsTrackRep Class Referenceabstract

Abstract base class for a track representation. More...

#include <AbsTrackRep.h>

Inheritance diagram for genfit::AbsTrackRep:
Collaboration diagram for genfit::AbsTrackRep:

Public Member Functions

 AbsTrackRep ()
 
 AbsTrackRep (int pdgCode, char propDir=0)
 
virtual ~AbsTrackRep ()
 
virtual AbsTrackRepclone () 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 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::MatStepgetSteps () 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 getTOF () const =0
 Get the time of flight [ns] of the last 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.
 
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!
 
AbsTrackRepoperator= (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_
 

Detailed Description

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.

Constructor & Destructor Documentation

◆ AbsTrackRep() [1/3]

genfit::AbsTrackRep::AbsTrackRep ( )

Definition at line 31 of file AbsTrackRep.cc.

31 :
32 pdgCode_(0), propDir_(0), debugLvl_(0)
33{
34 ;
35}
int pdgCode_
Particle code.
unsigned int debugLvl_
char propDir_
propagation direction (-1, 0, 1) -> (backward, auto, forward)

◆ AbsTrackRep() [2/3]

genfit::AbsTrackRep::AbsTrackRep ( int  pdgCode,
char  propDir = 0 
)

Definition at line 37 of file AbsTrackRep.cc.

37 :
38 pdgCode_(pdgCode), propDir_(propDir), debugLvl_(0)
39{
40 ;
41}

◆ ~AbsTrackRep()

virtual genfit::AbsTrackRep::~AbsTrackRep ( )
inlinevirtual

Definition at line 73 of file AbsTrackRep.h.

73{;}

◆ AbsTrackRep() [3/3]

genfit::AbsTrackRep::AbsTrackRep ( const AbsTrackRep rep)
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.

43 :
44 TObject(rep), pdgCode_(rep.pdgCode_), propDir_(rep.propDir_), debugLvl_(rep.debugLvl_)
45{
46 ;
47}

Member Function Documentation

◆ calcJacobianNumerically()

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.

107 {
108
109 // Find the transport matrix for track propagation from origState to destPlane
110 // I.e. this finds
111 // d stateDestPlane / d origState |_origState
112
113 jacobian.ResizeTo(getDim(), getDim());
114
115 // no science behind these values, I verified that forward and
116 // backward propagation yield inverse matrices to good
117 // approximation. In order to avoid bad roundoff errors, the actual
118 // step taken is determined below, separately for each direction.
119 const double defaultStepX = 1.E-5;
120 double stepX;
121
122 // Calculate derivative for all three dimensions successively.
123 // The algorithm follows the one in TF1::Derivative() :
124 // df(x) = (4 D(h/2) - D(h)) / 3
125 // with D(h) = (f(x + h) - f(x - h)) / (2 h).
126 //
127 // Could perhaps do better by also using f(x) which would be stB.
128 TVectorD rightShort(getDim()), rightFull(getDim());
129 TVectorD leftShort(getDim()), leftFull(getDim());
130 for (size_t i = 0; i < getDim(); ++i) {
131 {
132 genfit::StateOnPlane stateCopy(origState);
133 double temp = stateCopy.getState()(i) + defaultStepX / 2;
134 // Find the actual size of the step, which will differ from
135 // defaultStepX due to roundoff. This is the step-size we will
136 // use for this direction. Idea taken from Numerical Recipes,
137 // 3rd ed., section 5.7.
138 //
139 // Note that if a number is exactly representable, it's double
140 // will also be exact. Outside denormals, this also holds for
141 // halving. Unless the exponent changes (which it only will in
142 // the vicinity of zero) adding or subtracing doesn't make a
143 // difference.
144 //
145 // We determine the roundoff error for the half-step. If this
146 // is exactly representable, the full step will also be.
147 stepX = 2 * (temp - stateCopy.getState()(i));
148 (stateCopy.getState())(i) = temp;
149 extrapolateToPlane(stateCopy, destPlane);
150 rightShort = stateCopy.getState();
151 }
152 {
153 genfit::StateOnPlane stateCopy(origState);
154 (stateCopy.getState())(i) -= stepX / 2;
155 extrapolateToPlane(stateCopy, destPlane);
156 leftShort = stateCopy.getState();
157 }
158 {
159 genfit::StateOnPlane stateCopy(origState);
160 (stateCopy.getState())(i) += stepX;
161 extrapolateToPlane(stateCopy, destPlane);
162 rightFull = stateCopy.getState();
163 }
164 {
165 genfit::StateOnPlane stateCopy(origState);
166 (stateCopy.getState())(i) -= stepX;
167 extrapolateToPlane(stateCopy, destPlane);
168 leftFull = stateCopy.getState();
169 }
170
171 // Calculate the derivatives for the individual components of
172 // the track parameters.
173 for (size_t j = 0; j < getDim(); ++j) {
174 double derivFull = (rightFull(j) - leftFull(j)) / 2 / stepX;
175 double derivShort = (rightShort(j) - leftShort(j)) / stepX;
176
177 jacobian(j, i) = 1./3.*(4*derivShort - derivFull);
178 }
179 }
180}
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,...
virtual unsigned int getDim() const =0
Get the dimension of the state vector used by the track representation.
A state with arbitrary dimension defined in a DetPlane.
int i
Definition ShipAna.py:86

◆ clone()

virtual AbsTrackRep * genfit::AbsTrackRep::clone ( ) const
pure virtual

Clone the trackRep.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ extrapolateBy()

virtual double genfit::AbsTrackRep::extrapolateBy ( StateOnPlane state,
double  step,
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
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, and genfit::RKTrackRep.

◆ extrapolateToCylinder()

virtual double genfit::AbsTrackRep::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
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, and genfit::RKTrackRep.

◆ extrapolateToLine() [1/2]

virtual double genfit::AbsTrackRep::extrapolateToLine ( StateOnPlane state,
const TVector3 &  linePoint,
const TVector3 &  lineDirection,
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
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, genfit::RKTrackRep, genfit::RKTrackRep, and genfit::RKTrackRep.

◆ extrapolateToLine() [2/2]

virtual double genfit::AbsTrackRep::extrapolateToLine ( StateOnPlane state,
const TVector3 &  point1,
const TVector3 &  point2,
TVector3 &  poca,
TVector3 &  dirInPoca,
TVector3 &  poca_onwire,
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
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.

Reimplemented in genfit::RKTrackRep, and genfit::RKTrackRep.

Definition at line 120 of file AbsTrackRep.h.

127 {
128 TVector3 wireDir(point2 - point1);
129 wireDir.Unit();
130 double retval = this->extrapolateToLine(state, point1, wireDir, stopAtBoundary, calcJacobianNoise);
131 poca = this->getPos(state);
132 dirInPoca = this->getMom(state);
133 dirInPoca.Unit();
134
135 poca_onwire = point1 + wireDir*((poca - point1)*wireDir);
136
137 return retval;
138 }
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,...
virtual TVector3 getMom(const StateOnPlane &state) const =0
Get the cartesian momentum vector of a state.
virtual TVector3 getPos(const StateOnPlane &state) const =0
Get the cartesian position of a state.

◆ extrapolateToMeasurement()

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.

53 {
54
55 return this->extrapolateToPlane(state, measurement->constructPlane(state), stopAtBoundary, calcJacobianNoise);
56}

◆ extrapolateToPlane()

virtual double genfit::AbsTrackRep::extrapolateToPlane ( StateOnPlane state,
const genfit::SharedPlanePtr plane,
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
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, and genfit::RKTrackRep.

◆ extrapolateToPoint() [1/2]

virtual double genfit::AbsTrackRep::extrapolateToPoint ( StateOnPlane state,
const TVector3 &  point,
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
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, and genfit::RKTrackRep.

◆ extrapolateToPoint() [2/2]

virtual double genfit::AbsTrackRep::extrapolateToPoint ( StateOnPlane state,
const TVector3 &  point,
const TMatrixDSym &  G,
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
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, and genfit::RKTrackRep.

◆ extrapolateToSphere()

virtual double genfit::AbsTrackRep::extrapolateToSphere ( StateOnPlane state,
double  radius,
const TVector3 &  point = TVector3(0., 0., 0.),
bool  stopAtBoundary = false,
bool  calcJacobianNoise = false 
) const
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, and genfit::RKTrackRep.

◆ get6DCov()

virtual TMatrixDSym genfit::AbsTrackRep::get6DCov ( const MeasuredStateOnPlane state) const
pure virtual

Get the 6D covariance.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ get6DState()

TVectorD genfit::AbsTrackRep::get6DState ( const StateOnPlane state) const
virtual

Get the 6D state vector (x, y, z, p_x, p_y, p_z).

Definition at line 59 of file AbsTrackRep.cc.

59 {
60 TVector3 pos, mom;
61 getPosMom(state, pos, mom);
62
63 TVectorD stateVec(6);
64
65 stateVec(0) = pos.X();
66 stateVec(1) = pos.Y();
67 stateVec(2) = pos.Z();
68
69 stateVec(3) = mom.X();
70 stateVec(4) = mom.Y();
71 stateVec(5) = mom.Z();
72
73 return stateVec;
74}
virtual void getPosMom(const StateOnPlane &state, TVector3 &pos, TVector3 &mom) const =0
Get cartesian position and momentum vector of a state.

◆ get6DStateCov()

void genfit::AbsTrackRep::get6DStateCov ( const MeasuredStateOnPlane state,
TVectorD &  stateVec,
TMatrixDSym &  cov 
) const
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.

77 {
78 TVector3 pos, mom;
79 getPosMomCov(state, pos, mom, cov);
80
81 stateVec.ResizeTo(6);
82
83 stateVec(0) = pos.X();
84 stateVec(1) = pos.Y();
85 stateVec(2) = pos.Z();
86
87 stateVec(3) = mom.X();
88 stateVec(4) = mom.Y();
89 stateVec(5) = mom.Z();
90}
virtual void getPosMomCov(const MeasuredStateOnPlane &state, TVector3 &pos, TVector3 &mom, TMatrixDSym &cov) const =0
Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance.

◆ getBackwardJacobianAndNoise()

virtual void genfit::AbsTrackRep::getBackwardJacobianAndNoise ( TMatrixD &  jacobian,
TMatrixDSym &  noise,
TVectorD &  deltaState 
) const
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, and genfit::RKTrackRep.

◆ getCharge()

virtual double genfit::AbsTrackRep::getCharge ( const StateOnPlane state) const
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::RKTrackRep.

◆ getDim()

virtual unsigned int genfit::AbsTrackRep::getDim ( ) const
pure virtual

Get the dimension of the state vector used by the track representation.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ getDir()

TVector3 genfit::AbsTrackRep::getDir ( const StateOnPlane state) const
inline

Get the direction vector of a state.

Definition at line 230 of file AbsTrackRep.h.

230{return getMom(state).Unit();}

◆ getForwardJacobianAndNoise()

virtual void genfit::AbsTrackRep::getForwardJacobianAndNoise ( TMatrixD &  jacobian,
TMatrixDSym &  noise,
TVectorD &  deltaState 
) const
pure virtual

Get the jacobian and noise matrix of the last extrapolation.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ getMass()

double genfit::AbsTrackRep::getMass ( const StateOnPlane state) const

Get tha particle mass in GeV/c^2.

Definition at line 100 of file AbsTrackRep.cc.

100 {
101 return TDatabasePDG::Instance()->GetParticle(pdgCode_)->Mass();
102}

◆ getMom()

virtual TVector3 genfit::AbsTrackRep::getMom ( const StateOnPlane state) const
pure virtual

Get the cartesian momentum vector of a state.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ getMomMag()

virtual double genfit::AbsTrackRep::getMomMag ( const StateOnPlane state) const
pure virtual

get the magnitude of the momentum in GeV.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ getMomVar()

virtual double genfit::AbsTrackRep::getMomVar ( const MeasuredStateOnPlane state) const
pure virtual

get the variance of the absolute value of the momentum .

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ getPDG()

int genfit::AbsTrackRep::getPDG ( ) const
inline

Get the pdg code.

Definition at line 256 of file AbsTrackRep.h.

256{return pdgCode_;}

◆ getPDGCharge()

double genfit::AbsTrackRep::getPDGCharge ( ) const

Get the charge of the particle of the pdg code.

Definition at line 93 of file AbsTrackRep.cc.

93 {
94 TParticlePDG* particle = TDatabasePDG::Instance()->GetParticle(pdgCode_);
95 assert(particle != NULL);
96 return particle->Charge()/(3.);
97}

◆ getPos()

virtual TVector3 genfit::AbsTrackRep::getPos ( const StateOnPlane state) const
pure virtual

Get the cartesian position of a state.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ getPosDir()

void genfit::AbsTrackRep::getPosDir ( const StateOnPlane state,
TVector3 &  pos,
TVector3 &  dir 
) const
inline

Get cartesian position and direction vector of a state.

Definition at line 236 of file AbsTrackRep.h.

236{getPosMom(state, pos, dir); dir.SetMag(1.);}

◆ getPosMom()

virtual void genfit::AbsTrackRep::getPosMom ( const StateOnPlane state,
TVector3 &  pos,
TVector3 &  mom 
) const
pure virtual

Get cartesian position and momentum vector of a state.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ getPosMomCov()

virtual void genfit::AbsTrackRep::getPosMomCov ( const MeasuredStateOnPlane state,
TVector3 &  pos,
TVector3 &  mom,
TMatrixDSym &  cov 
) const
pure virtual

Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ getPropDir()

char genfit::AbsTrackRep::getPropDir ( ) const
inline

Get propagation direction. (-1, 0, 1) -> (backward, auto, forward).

Definition at line 272 of file AbsTrackRep.h.

272{return propDir_;}

◆ getQop()

virtual double genfit::AbsTrackRep::getQop ( const StateOnPlane state) const
pure virtual

Get charge over momentum.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ getRadiationLenght()

virtual double genfit::AbsTrackRep::getRadiationLenght ( ) const
pure virtual

Get the accumulated X/X0 (path / radiation length) of the material crossed in the last extrapolation.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ getSteps()

virtual std::vector< genfit::MatStep > genfit::AbsTrackRep::getSteps ( ) const
pure virtual

Get stepsizes and material properties of crossed materials of the last extrapolation.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ getTOF()

virtual double genfit::AbsTrackRep::getTOF ( ) const
pure virtual

Get the time of flight [ns] of the last extrapolation.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ isSame()

virtual bool genfit::AbsTrackRep::isSame ( const AbsTrackRep other)
pure virtual

check if other is of same type (e.g. RKTrackRep) and has same pdg code.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ isSameType()

virtual bool genfit::AbsTrackRep::isSameType ( const AbsTrackRep other)
pure virtual

check if other is of same type (e.g. RKTrackRep).

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ operator=()

AbsTrackRep & genfit::AbsTrackRep::operator= ( const AbsTrackRep )
protected

protect from calling assignment operator from outside the class. Use clone() instead!

◆ Print()

void genfit::AbsTrackRep::Print ( const Option_t *  = "") const
virtual

Definition at line 194 of file AbsTrackRep.cc.

194 {
195 std::cout << "genfit::AbsTrackRep, pdgCode = " << pdgCode_ << ". PropDir = " << (int)propDir_ << "\n";
196}

◆ setChargeSign()

virtual void genfit::AbsTrackRep::setChargeSign ( StateOnPlane state,
double  charge 
) const
pure virtual

Set the sign of the charge according to charge.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ setDebugLvl()

virtual void genfit::AbsTrackRep::setDebugLvl ( unsigned int  lvl = 1)
inlinevirtual

Definition at line 332 of file AbsTrackRep.h.

332{debugLvl_ = lvl;}

◆ setPosMom() [1/2]

virtual void genfit::AbsTrackRep::setPosMom ( StateOnPlane state,
const TVector3 &  pos,
const TVector3 &  mom 
) const
pure virtual

Set position and momentum of state.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ setPosMom() [2/2]

virtual void genfit::AbsTrackRep::setPosMom ( StateOnPlane state,
const TVectorD &  state6 
) const
pure virtual

Set position and momentum of state.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ setPosMomCov() [1/2]

virtual void genfit::AbsTrackRep::setPosMomCov ( MeasuredStateOnPlane state,
const TVector3 &  pos,
const TVector3 &  mom,
const TMatrixDSym &  cov6x6 
) const
pure virtual

Set position, momentum and covariance of state.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ setPosMomCov() [2/2]

virtual void genfit::AbsTrackRep::setPosMomCov ( MeasuredStateOnPlane state,
const TVectorD &  state6,
const TMatrixDSym &  cov6x6 
) const
pure virtual

Set position, momentum and covariance of state.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ setPosMomErr()

virtual void genfit::AbsTrackRep::setPosMomErr ( MeasuredStateOnPlane state,
const TVector3 &  pos,
const TVector3 &  mom,
const TVector3 &  posErr,
const TVector3 &  momErr 
) const
pure virtual

Set position and momentum and error of state.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ setPropDir()

void genfit::AbsTrackRep::setPropDir ( int  dir)
inline

Set propagation direction. (-1, 0, 1) -> (backward, auto, forward).

Definition at line 317 of file AbsTrackRep.h.

317 {
318 if (dir>0) propDir_ = 1;
319 else if (dir<0) propDir_ = -1;
320 else propDir_ = 0;
321 };

◆ setQop()

virtual void genfit::AbsTrackRep::setQop ( StateOnPlane state,
double  qop 
) const
pure virtual

Set charge/momentum.

Implemented in genfit::RKTrackRep, and genfit::RKTrackRep.

◆ switchPDGSign()

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.

183 {
184 TParticlePDG* particle = TDatabasePDG::Instance()->GetParticle(-pdgCode_);
185 if(particle != NULL) {
186 pdgCode_ *= -1;
187 return true;
188 }
189 return false;
190}

◆ switchPropDir()

void genfit::AbsTrackRep::switchPropDir ( )
inline

Switch propagation direction. Has no effect if propDir_ is set to 0.

Definition at line 324 of file AbsTrackRep.h.

324{propDir_ = -1*propDir_;}

Member Data Documentation

◆ debugLvl_

unsigned int genfit::AbsTrackRep::debugLvl_
protected

Definition at line 349 of file AbsTrackRep.h.

◆ pdgCode_

int genfit::AbsTrackRep::pdgCode_
protected

Particle code.

Definition at line 345 of file AbsTrackRep.h.

◆ propDir_

char genfit::AbsTrackRep::propDir_
protected

propagation direction (-1, 0, 1) -> (backward, auto, forward)

Definition at line 347 of file AbsTrackRep.h.


The documentation for this class was generated from the following files: