Belle II Software development
CDCTrajectory3D Class Reference

Particle full three dimensional trajectory. More...

#include <CDCTrajectory3D.h>

Public Member Functions

 CDCTrajectory3D ()
 Default constructor for ROOT compatibility.
 
 CDCTrajectory3D (const UncertainHelix &helix)
 Constructs a trajectory from a helix with reference point equivalent to the origin.
 
 CDCTrajectory3D (const Helix &helix)
 conversion constructor to make that one stupid test work
 
 CDCTrajectory3D (const Vector3D &localOrigin, const UncertainHelix &localHelix, double flightTime=NAN)
 Constructs a trajectory from a local helix taken as relative to the given origin.
 
 CDCTrajectory3D (const CDCTrajectory2D &trajectory2D, const CDCTrajectorySZ &trajectorySZ)
 Construct a three dimensional trajectory from a two dimensional circular trajectory and sz linear trajectory.
 
 CDCTrajectory3D (const CDCTrajectory2D &trajectory2D)
 Construct a trajectory from a two dimensional circular trajectory filling the remaining two parameters and covariance matrix with default values.
 
 CDCTrajectory3D (const Vector3D &pos3D, double time, const Vector3D &mom3D, double charge, double bZ)
 Construct a trajectory with given start point, momentum at the start point and given charge.
 
 CDCTrajectory3D (const Vector3D &pos3D, double time, const Vector3D &mom3D, double charge)
 Construct a trajectory with given start point, momentum at the start point and given charge.
 
 CDCTrajectory3D (const MCParticle &mcParticle, double bZ)
 Construct a trajectory from the MCParticles vertex and momentum.
 
 CDCTrajectory3D (const MCParticle &mcParticle)
 Construct a trajectory from the MCParticles vertex and momentum.
 
 CDCTrajectory3D (const genfit::TrackCand &gfTrackCand, double bZ)
 Construct a trajectory by extracting the seed position of the genfit::TrackCand.
 
 CDCTrajectory3D (const genfit::TrackCand &gfTrackCand)
 Construct a trajectory by extracting the seed position of the genfit::TrackCand.
 
bool isInvalid () const
 Checks if the trajectory is already set to a valid value.
 
bool isFitted () const
 Checks if the trajectory has already been set to a valid value.
 
void clear ()
 Clears all information from this trajectoy.
 
bool fillInto (genfit::TrackCand &trackCand) const
 Copies the trajectory information to the Genfit track candidate.
 
bool fillInto (genfit::TrackCand &gfTrackCand, double bZ) const
 Copies the trajectory information to the Genfit track candidate.
 
CovarianceMatrix< 6 > getCartesianCovariance (double bZ) const
 Convert the helix parameters to the cartesian coordinates x,y,z,px,py,pz.
 
void reverse ()
 Reverses the trajectory in place.
 
CDCTrajectory3D reversed () const
 Returns the reverse trajectory as a copy.
 
Vector3D reconstruct3D (const WireLine &wireLine, double distance=0.0) const
 Gives the three dimensional point which is on the dirft circle away from the wire line.
 
double calcArcLength2D (const Vector3D &point) const
 Calculate the travel distance from the start position of the trajectory.
 
double getArcLength2DPeriod () const
 Getter for the arc length for one round trip around the trajectory.
 
double shiftPeriod (int nPeriods)
 Adjusts the z0 to the one that lies n periods forward.
 
ESign getChargeSign () const
 Gets the charge sign of the trajectory.
 
double getAbsMom3D (double bZ) const
 Get the estimation for the absolute value of the transvers momentum.
 
double getAbsMom3D () const
 Get the estimation for the absolute value of the transvers momentum.
 
Vector3D getMom3DAtSupport (const double bZ) const
 Get the momentum at the start point of the trajectory.
 
Vector3D getMom3DAtSupport () const
 Get the momentum at the start point of the trajectory.
 
Vector3D getFlightDirection3DAtSupport () const
 Get the unit momentum at the start point of the trajectory.
 
Vector3D getSupport () const
 Getter for the support point of the trajectory in global coordinates, where arcLength2D = 0.
 
Vector3D getGlobalPerigee () const
 Getter for the closest approach on the trajectory to the global origin.
 
Vector2D getGlobalCenter () const
 Getter for the center of the helix in global coordinates.
 
bool isCurler (double factor=1) const
 Checks if the trajectory leaves the outer radius of the CDC times the given tolerance factor.
 
double getMaximalCylindricalR () const
 Getter for the maximal distance from the origin.
 
double getMinimalCylindricalR () const
 Getter for the minimal distance from the origin.
 
double getGlobalImpact () const
 Getter for the signed impact parameter of the trajectory.
 
CDCTrajectory2D getTrajectory2D () const
 Getter for the two dimensional trajectory.
 
CDCTrajectorySZ getTrajectorySZ () const
 Getter for the sz trajectory.
 
PerigeeCircle getGlobalCircle () const
 Getter for the circle in global coordinates.
 
UncertainPerigeeCircle getLocalCircle () const
 Getter for the circle in local coordinates.
 
UncertainSZLine getLocalSZLine () const
 Getter for the sz line starting from the local origin.
 
double getLocalCovariance (EHelixParameter iRow, EHelixParameter iCol) const
 Getter for an individual element of the covariance matrix of the local helix parameters.
 
double getLocalVariance (EHelixParameter i) const
 Getter for an individual diagonal element of the covariance matrix of the local helix parameters.
 
double getTanLambda () const
 Getter for the slope of z over the transverse travel distance s.
 
double getCurvatureXY () const
 Getter for the curvature as seen from the xy projection.
 
double getPValue () const
 Getter for p-value.
 
double getChi2 () const
 Getter for the chi2 value of the fit.
 
void setChi2 (const double chi2)
 Setter for the chi square value of the helix fit.
 
size_t getNDF () const
 Getter for the number of degrees of freedom of the helix fit.
 
void setNDF (std::size_t ndf)
 Setter for the number of degrees of freedom of the helix fit.
 
const UncertainHelixgetLocalHelix () const
 Getter for the helix in local coordinates.
 
void setLocalHelix (const UncertainHelix &localHelix)
 Setter for the helix that describes the trajectory in local coordinates.
 
const Vector3DgetLocalOrigin () const
 Getter for the origin of the local coordinate system.
 
double setLocalOrigin (const Vector3D &localOrigin)
 Setter for the origin of the local coordinate system.
 
double getFlightTime () const
 Getter for the time when the particle reached the support point position.
 
void setFlightTime (double flightTime)
 Setter for the time when the particle reached the support point position.
 

Private Attributes

Vector3D m_localOrigin
 Memory for local coordinate origin of the circle representing the trajectory in global coordinates.
 
UncertainHelix m_localHelix
 Memory for the generalized circle describing the trajectory in coordinates from the local origin.
 
double m_flightTime = NAN
 Memory for the estimation of the time at which the particle arrived at the support point.
 

Detailed Description

Particle full three dimensional trajectory.

Definition at line 45 of file CDCTrajectory3D.h.

Constructor & Destructor Documentation

◆ CDCTrajectory3D() [1/12]

CDCTrajectory3D ( )
inline

Default constructor for ROOT compatibility.

Definition at line 49 of file CDCTrajectory3D.h.

51 , m_localHelix()
52 {
53 }
Vector3D m_localOrigin
Memory for local coordinate origin of the circle representing the trajectory in global coordinates.
UncertainHelix m_localHelix
Memory for the generalized circle describing the trajectory in coordinates from the local origin.

◆ CDCTrajectory3D() [2/12]

CDCTrajectory3D ( const UncertainHelix helix)
inlineexplicit

Constructs a trajectory from a helix with reference point equivalent to the origin.

Definition at line 56 of file CDCTrajectory3D.h.

57 : m_localOrigin(0.0, 0.0, 0.0)
58 , m_localHelix(helix)
59 , m_flightTime(0.0)
60 {
61 }
double m_flightTime
Memory for the estimation of the time at which the particle arrived at the support point.

◆ CDCTrajectory3D() [3/12]

CDCTrajectory3D ( const Helix helix)
inlineexplicit

conversion constructor to make that one stupid test work

Definition at line 64 of file CDCTrajectory3D.h.

64 : CDCTrajectory3D(UncertainHelix(helix))
65 {
66 }
CDCTrajectory3D()
Default constructor for ROOT compatibility.

◆ CDCTrajectory3D() [4/12]

CDCTrajectory3D ( const Vector3D localOrigin,
const UncertainHelix localHelix,
double  flightTime = NAN 
)
inline

Constructs a trajectory from a local helix taken as relative to the given origin.

Definition at line 69 of file CDCTrajectory3D.h.

72 : m_localOrigin(localOrigin)
73 , m_localHelix(localHelix)
74 , m_flightTime(flightTime)
75 {
76 }

◆ CDCTrajectory3D() [5/12]

CDCTrajectory3D ( const CDCTrajectory2D trajectory2D,
const CDCTrajectorySZ trajectorySZ 
)

Construct a three dimensional trajectory from a two dimensional circular trajectory and sz linear trajectory.

Definition at line 47 of file CDCTrajectory3D.cc.

49 : m_localOrigin(trajectory2D.getLocalOrigin())
50 , m_localHelix(trajectory2D.getLocalCircle(), trajectorySZ.getSZLine())
51 , m_flightTime(trajectory2D.getFlightTime())
52{
53}
double getFlightTime() const
Getter for the time when the particle reached the support point position.
const Vector2D & getLocalOrigin() const
Getter for the origin of the local coordinate system.
const UncertainPerigeeCircle & getLocalCircle() const
Getter for the cirlce in local coordinates.
const UncertainSZLine & getSZLine() const
Getter for the line in sz space.

◆ CDCTrajectory3D() [6/12]

CDCTrajectory3D ( const CDCTrajectory2D trajectory2D)
explicit

Construct a trajectory from a two dimensional circular trajectory filling the remaining two parameters and covariance matrix with default values.

Definition at line 55 of file CDCTrajectory3D.cc.

57{
58}
static CDCTrajectorySZ basicAssumption()
Constucts a basic assumption, what the z0 start position and the sz slope are, including some broad v...

◆ CDCTrajectory3D() [7/12]

CDCTrajectory3D ( const Vector3D pos3D,
double  time,
const Vector3D mom3D,
double  charge,
double  bZ 
)

Construct a trajectory with given start point, momentum at the start point and given charge.

Additionally this can takes an explicit bZ value instead of a field value from the instance BFieldMap.

Definition at line 60 of file CDCTrajectory3D.cc.

65 : m_localOrigin(pos3D)
67 mom3D.xy().unit(),
68 0.0,
69 mom3D.cotTheta(),
70 0.0)
71 , m_flightTime(time)
72{
73}
static double absMom2DToCurvature(double absMom2D, double charge, double bZ)
Conversion helper for momenta to two dimensional curvature.
Vector2D unit() const
Returns a unit vector colaligned with this.
Definition: Vector2D.h:321
double norm() const
Calculates the length of the vector.
Definition: Vector2D.h:175
const Vector2D & xy() const
Getter for the xy projected vector ( reference ! )
Definition: Vector3D.h:508
double cotTheta() const
Getter for the cotangent of the polar angle.
Definition: Vector3D.h:558

◆ CDCTrajectory3D() [8/12]

CDCTrajectory3D ( const Vector3D pos3D,
double  time,
const Vector3D mom3D,
double  charge 
)

Construct a trajectory with given start point, momentum at the start point and given charge.

Definition at line 75 of file CDCTrajectory3D.cc.

79 : CDCTrajectory3D(pos3D, time, mom3D, charge, CDCBFieldUtil::getBFieldZ(pos3D))
80{
81}
static double getBFieldZ()
Getter for the signed magnetic field stength in z direction at the origin ( in Tesla )

◆ CDCTrajectory3D() [9/12]

CDCTrajectory3D ( const MCParticle mcParticle,
double  bZ 
)

Construct a trajectory from the MCParticles vertex and momentum.

Definition at line 83 of file CDCTrajectory3D.cc.

85 mcParticle.getProductionTime(),
86 Vector3D{mcParticle.getMomentum()},
87 mcParticle.getCharge(),
88 bZ)
89{
90}
ROOT::Math::XYZVector getProductionVertex() const
Return production vertex position.
Definition: MCParticle.h:189
float getCharge() const
Return the particle charge defined in TDatabasePDG.
Definition: MCParticle.cc:36
float getProductionTime() const
Return production time in ns.
Definition: MCParticle.h:159
ROOT::Math::XYZVector getMomentum() const
Return momentum.
Definition: MCParticle.h:198
A three dimensional vector.
Definition: Vector3D.h:33

◆ CDCTrajectory3D() [10/12]

CDCTrajectory3D ( const MCParticle mcParticle)
explicit

Construct a trajectory from the MCParticles vertex and momentum.

Definition at line 92 of file CDCTrajectory3D.cc.

94 mcParticle.getProductionTime(),
95 Vector3D{mcParticle.getMomentum()},
96 mcParticle.getCharge())
97{
98}

◆ CDCTrajectory3D() [11/12]

CDCTrajectory3D ( const genfit::TrackCand &  gfTrackCand,
double  bZ 
)

Construct a trajectory by extracting the seed position of the genfit::TrackCand.

Definition at line 100 of file CDCTrajectory3D.cc.

101 : CDCTrajectory3D(Vector3D{gfTrackCand.getPosSeed()},
102 gfTrackCand.getTimeSeed(),
103 Vector3D{gfTrackCand.getMomSeed()},
104 gfTrackCand.getChargeSeed(),
105 bZ)
106{
107 // Maybe push these out of this function:
108 // Indices of the cartesian coordinates
109 const int iX = 0;
110 const int iY = 1;
111 const int iZ = 2;
112 const int iPx = 3;
113 const int iPy = 4;
114 const int iPz = 5;
115
116 const TMatrixDSym& seedCov = gfTrackCand.getCovSeed();
117 CovarianceMatrix<6> cov6 = TMatrixConversion::fromTMatrix<6>(seedCov);
118
119 // 1. Rotate to a system where phi0 = 0
120 JacobianMatrix<6, 6> jacobianRot = JacobianMatrixUtil::zero<6, 6>();
121
122 const double px = gfTrackCand.getStateSeed()[iPx];
123 const double py = gfTrackCand.getStateSeed()[iPy];
124 const double pt = hypot2(px, py);
125
126 const double cosPhi0 = px / pt;
127 const double sinPhi0 = py / pt;
128
129 // Passive rotation matrix by phi0:
130 jacobianRot(iX, iX) = cosPhi0;
131 jacobianRot(iX, iY) = sinPhi0;
132 jacobianRot(iY, iX) = -sinPhi0;
133 jacobianRot(iY, iY) = cosPhi0;
134 jacobianRot(iZ, iZ) = 1.0;
135
136 jacobianRot(iPx, iPx) = cosPhi0;
137 jacobianRot(iPx, iPy) = sinPhi0;
138 jacobianRot(iPy, iPx) = -sinPhi0;
139 jacobianRot(iPy, iPy) = cosPhi0;
140 jacobianRot(iPz, iPz) = 1.0;
141
142 // Apply transformation inplace
143 CovarianceMatrixUtil::transport(jacobianRot, cov6);
144
145 // 2. Translate to perigee parameters
146 JacobianMatrix<5, 6> jacobianReduce = JacobianMatrixUtil::zero<5, 6>();
147
148 const double invPt = 1 / pt;
149 const double invPtSquared = invPt * invPt;
150 const double pz = gfTrackCand.getStateSeed()[iPz];
151 const double alpha = CDCBFieldUtil::getAlphaFromBField(bZ);
152 const double charge = gfTrackCand.getChargeSeed();
153
154 using namespace NHelixParameterIndices;
155 jacobianReduce(c_Curv, iPx) = charge * invPtSquared / alpha ;
156 jacobianReduce(c_Phi0, iPy) = invPt;
157 jacobianReduce(c_I, iY) = 1;
158 jacobianReduce(c_TanL, iPx) = - pz * invPtSquared;
159 jacobianReduce(c_TanL, iPz) = invPt;
160 jacobianReduce(c_Z0, iZ) = 1;
161 // Note the column corresponding to iX is completely zero as expectable.
162
163 CovarianceMatrix<5> cov5 = CovarianceMatrixUtil::transported(jacobianReduce, cov6);
164 const HelixCovariance& helixCovariance = cov5;
165
166 // The covariance should now be the correct 5x5 covariance matrix.
167 m_localHelix.setHelixCovariance(helixCovariance);
168}
static double getAlphaFromBField(double bField)
Translater from magnetic field strength in Tesla to the alpha value.
A matrix implementation to be used as an interface typ through out the track finder.
Definition: PlainMatrix.h:40
void setHelixCovariance(const HelixCovariance &helixCovariance)
Setter for the whole covariance matrix of the perigee parameters.
double charge(int pdgCode)
Returns electric charge of a particle with given pdg code.
Definition: EvtPDLUtil.cc:44
static CovarianceMatrix< M > transported(const JacobianMatrix< M, N > &jacobian, const CovarianceMatrix< N > &cov)
Return a copy of the covariance matrix transported with the given jacobian matrix.
static void transport(const JacobianMatrix< N, N > &jacobian, CovarianceMatrix< N > &cov)
Transport the covariance matrix inplace with the given jacobian matrix.

◆ CDCTrajectory3D() [12/12]

CDCTrajectory3D ( const genfit::TrackCand &  gfTrackCand)
explicit

Construct a trajectory by extracting the seed position of the genfit::TrackCand.

Definition at line 170 of file CDCTrajectory3D.cc.

171 : CDCTrajectory3D(gfTrackCand, CDCBFieldUtil::getBFieldZ(Vector3D{gfTrackCand.getPosSeed()}))
172{
173}

Member Function Documentation

◆ calcArcLength2D()

double calcArcLength2D ( const Vector3D point) const
inline

Calculate the travel distance from the start position of the trajectory.

Returns the travel distance on the trajectory from the start point to
the given point. This is subjected to a discontinuity at the far point
of the circle. Hence the value return is in the range from -pi*radius to pi*radius
If you have a heavily curling track you have care about the feasibility of this
calculation.

Definition at line 178 of file CDCTrajectory3D.h.

179 {
180 return getLocalHelix()->circleXY().arcLengthTo((point - getLocalOrigin()).xy());
181 }
const UncertainHelix & getLocalHelix() const
Getter for the helix in local coordinates.
const Vector3D & getLocalOrigin() const
Getter for the origin of the local coordinate system.
const PerigeeCircle & circleXY() const
Getter for the projection into xy space.
Definition: Helix.h:320
double arcLengthTo(const Vector2D &point) const
Calculates the arc length between the perigee and the given point.

◆ clear()

void clear ( )
inline

Clears all information from this trajectoy.

Definition at line 130 of file CDCTrajectory3D.h.

131 {
132 m_localOrigin.set(0.0, 0.0, 0.0);
134 m_flightTime = NAN;
135 }
void invalidate()
Sets all circle parameters to zero and the covariance matrix to something noninformative.
void set(const double first, const double second, const double third)
Setter for all three coordinates.
Definition: Vector3D.h:520

◆ fillInto() [1/2]

bool fillInto ( genfit::TrackCand &  gfTrackCand,
double  bZ 
) const

Copies the trajectory information to the Genfit track candidate.

Definition at line 266 of file CDCTrajectory3D.cc.

267{
268 // Set the start parameters
269 Vector3D position = getSupport();
270 Vector3D momentum = bZ == 0 ? getFlightDirection3DAtSupport() : getMom3DAtSupport(bZ);
272
273 // Do not propagate invalid fits, signal that the fit is invalid to the caller.
274 if (not ESignUtil::isValid(charge) or momentum.hasNAN() or position.hasNAN()) {
275 return false;
276 }
277
278 gfTrackCand.setPosMomSeed(position, momentum, charge);
279
280 const CovarianceMatrix<6> cov6 = calculateCovarianceMatrix(getLocalHelix(), momentum, charge, bZ);
281 TMatrixDSym covSeed = TMatrixConversion::toTMatrix(cov6);
282
283 gfTrackCand.setCovSeed(covSeed);
284
285 return true;
286}
ESign getChargeSign() const
Gets the charge sign of the trajectory.
Vector3D getFlightDirection3DAtSupport() const
Get the unit momentum at the start point of the trajectory.
Vector3D getMom3DAtSupport() const
Get the momentum at the start point of the trajectory.
Vector3D getSupport() const
Getter for the support point of the trajectory in global coordinates, where arcLength2D = 0.
bool hasNAN() const
Checks if one of the coordinates is NAN.
Definition: Vector3D.h:165
ESign
Enumeration for the distinct sign values of floating point variables.
Definition: ESign.h:27
bool isValid(ESign s)
Returns true if sign is ESign::c_Plus, ESign::c_Minus or ESign::c_Zero.
Definition: ESign.h:46
static TMatrixDSym toTMatrix(const CovarianceMatrix< N > &cov)
Translate the covariance matrix to the TMatrix representation.

◆ fillInto() [2/2]

bool fillInto ( genfit::TrackCand &  trackCand) const

Copies the trajectory information to the Genfit track candidate.

Definition at line 260 of file CDCTrajectory3D.cc.

261{
262 Vector3D position = getSupport();
263 return fillInto(trackCand, CDCBFieldUtil::getBFieldZ(position));
264}
bool fillInto(genfit::TrackCand &trackCand) const
Copies the trajectory information to the Genfit track candidate.

◆ getAbsMom3D() [1/2]

double getAbsMom3D ( ) const

Get the estimation for the absolute value of the transvers momentum.

Definition at line 316 of file CDCTrajectory3D.cc.

317{
318 Vector3D position = getSupport();
319 double tanLambda = getLocalHelix()->tanLambda();
320 double factor2DTo3D = hypot2(1, tanLambda);
321 double curvatureXY = getLocalHelix()->curvatureXY();
322 double absMom2D = CDCBFieldUtil::curvatureToAbsMom2D(curvatureXY, position);
323 return factor2DTo3D * absMom2D;
324}
static double curvatureToAbsMom2D(double curvature, double bZ)
Conversion helper for two dimensional curvature to momenta.
double tanLambda() const
Getter for the proportinality factor from arc length in xy space to z.
Definition: Helix.h:240
double curvatureXY() const
Getter for the signed curvature in the xy projection.
Definition: Helix.h:204

◆ getAbsMom3D() [2/2]

double getAbsMom3D ( double  bZ) const

Get the estimation for the absolute value of the transvers momentum.

Definition at line 307 of file CDCTrajectory3D.cc.

308{
309 double tanLambda = getLocalHelix()->tanLambda();
310 double factor2DTo3D = hypot2(1, tanLambda);
311 double curvatureXY = getLocalHelix()->curvatureXY();
312 double absMom2D = CDCBFieldUtil::curvatureToAbsMom2D(curvatureXY, bZ);
313 return factor2DTo3D * absMom2D;
314}

◆ getArcLength2DPeriod()

double getArcLength2DPeriod ( ) const
inline

Getter for the arc length for one round trip around the trajectory.

Definition at line 184 of file CDCTrajectory3D.h.

185 {
187 }
double arcLength2DPeriod() const
Getter for the arc length of one trip around the helix.
Definition: Helix.h:265

◆ getCartesianCovariance()

CovarianceMatrix< 6 > getCartesianCovariance ( double  bZ) const

Convert the helix parameters to the cartesian coordinates x,y,z,px,py,pz.

Definition at line 288 of file CDCTrajectory3D.cc.

289{
290 // Set the start parameters
291 Vector3D momentum = bZ == 0 ? getFlightDirection3DAtSupport() : getMom3DAtSupport(bZ);
293 return calculateCovarianceMatrix(getLocalHelix(), momentum, charge, bZ);
294}

◆ getChargeSign()

ESign getChargeSign ( ) const

Gets the charge sign of the trajectory.

Definition at line 302 of file CDCTrajectory3D.cc.

303{
304 return CDCBFieldUtil::ccwInfoToChargeSign(getLocalHelix()->circleXY().orientation());
305}
static ESign ccwInfoToChargeSign(ERotation ccwInfo)
Conversion helper from clockwise or counterclockwise travel to the charge sign.

◆ getChi2()

double getChi2 ( ) const
inline

Getter for the chi2 value of the fit.

Definition at line 315 of file CDCTrajectory3D.h.

316 {
317 return getLocalHelix().chi2();
318 }
double chi2() const
Getter for the chi square value of the helix fit.

◆ getCurvatureXY()

double getCurvatureXY ( ) const
inline

Getter for the curvature as seen from the xy projection.

Definition at line 303 of file CDCTrajectory3D.h.

304 {
305 return getLocalHelix()->curvatureXY();
306 }

◆ getFlightDirection3DAtSupport()

Vector3D getFlightDirection3DAtSupport ( ) const
inline

Get the unit momentum at the start point of the trajectory.

Definition at line 220 of file CDCTrajectory3D.h.

221 {
222 return getLocalHelix()->tangential();
223 }
Vector3D tangential() const
Getter for the unit three dimensional tangential vector at the perigee point of the helix.
Definition: Helix.h:289

◆ getFlightTime()

double getFlightTime ( ) const
inline

Getter for the time when the particle reached the support point position.

Definition at line 364 of file CDCTrajectory3D.h.

365 {
366 return m_flightTime;
367 }

◆ getGlobalCenter()

Vector2D getGlobalCenter ( ) const
inline

Getter for the center of the helix in global coordinates.

Definition at line 239 of file CDCTrajectory3D.h.

240 {
241 return getLocalHelix()->centerXY() + m_localOrigin.xy();
242 }
Vector2D centerXY() const
Getter for the central point of the helix.
Definition: Helix.h:283

◆ getGlobalCircle()

PerigeeCircle getGlobalCircle ( ) const

Getter for the circle in global coordinates.

Definition at line 348 of file CDCTrajectory3D.cc.

349{
350 // Down cast since we do not necessarily wont the covariance matrix transformed as well
351 PerigeeCircle result(getLocalHelix()->circleXY());
352 result.passiveMoveBy(-getLocalOrigin().xy());
353 return result;
354}
Extension of the generalized circle also caching the perigee coordinates.
Definition: PerigeeCircle.h:36

◆ getGlobalImpact()

double getGlobalImpact ( ) const
inline

Getter for the signed impact parameter of the trajectory.

Definition at line 262 of file CDCTrajectory3D.h.

263 {
265 }
double distanceXY(const Vector2D &point) const
Calculates the distance of the line parallel to the z axes through the given point.
Definition: Helix.h:138

◆ getGlobalPerigee()

Vector3D getGlobalPerigee ( ) const
inline

Getter for the closest approach on the trajectory to the global origin.

Definition at line 233 of file CDCTrajectory3D.h.

234 {
236 }
Vector3D closestXY(const Vector2D &pointXY) const
Calculates the point on the helix with the smallest perpendicular (xy) distance.
Definition: Helix.h:125

◆ getLocalCircle()

UncertainPerigeeCircle getLocalCircle ( ) const

Getter for the circle in local coordinates.

Definition at line 356 of file CDCTrajectory3D.cc.

357{
359}
UncertainPerigeeCircle uncertainCircleXY() const
Projects the helix into the xy plain carrying over the relevant parts of the convariance matrix.

◆ getLocalCovariance()

double getLocalCovariance ( EHelixParameter  iRow,
EHelixParameter  iCol 
) const
inline

Getter for an individual element of the covariance matrix of the local helix parameters.

Definition at line 284 of file CDCTrajectory3D.h.

285 {
286 return getLocalHelix().covariance(iRow, iCol);
287 }
double covariance(const EHelixParameter &iRow, const EHelixParameter &iCol) const
Getter for individual elements of the covariance matrix.

◆ getLocalHelix()

const UncertainHelix & getLocalHelix ( ) const
inline

Getter for the helix in local coordinates.

Definition at line 339 of file CDCTrajectory3D.h.

340 {
341 return m_localHelix;
342 }

◆ getLocalOrigin()

const Vector3D & getLocalOrigin ( ) const
inline

Getter for the origin of the local coordinate system.

Definition at line 351 of file CDCTrajectory3D.h.

352 {
353 return m_localOrigin;
354 }

◆ getLocalSZLine()

UncertainSZLine getLocalSZLine ( ) const

Getter for the sz line starting from the local origin.

Definition at line 361 of file CDCTrajectory3D.cc.

362{
364}
UncertainSZLine uncertainSZLine() const
Reduces the helix to an sz line carrying over the relevant parts of the convariance matrix.

◆ getLocalVariance()

double getLocalVariance ( EHelixParameter  i) const
inline

Getter for an individual diagonal element of the covariance matrix of the local helix parameters.

Definition at line 291 of file CDCTrajectory3D.h.

292 {
293 return getLocalHelix().variance(i);
294 }
double variance(const EHelixParameter &i) const
Getter for individual diagonal elements of the covariance matrix.

◆ getMaximalCylindricalR()

double getMaximalCylindricalR ( ) const
inline

Getter for the maximal distance from the origin.

Definition at line 250 of file CDCTrajectory3D.h.

251 {
252 return std::fabs(getGlobalImpact() + 2 * getLocalHelix()->radiusXY());
253 }
double getGlobalImpact() const
Getter for the signed impact parameter of the trajectory.

◆ getMinimalCylindricalR()

double getMinimalCylindricalR ( ) const
inline

Getter for the minimal distance from the origin.

Definition at line 256 of file CDCTrajectory3D.h.

257 {
258 return std::fabs(getGlobalImpact());
259 }

◆ getMom3DAtSupport() [1/2]

Vector3D getMom3DAtSupport ( ) const
inline

Get the momentum at the start point of the trajectory.

Definition at line 214 of file CDCTrajectory3D.h.

215 {
217 }
double getAbsMom3D() const
Get the estimation for the absolute value of the transvers momentum.

◆ getMom3DAtSupport() [2/2]

Vector3D getMom3DAtSupport ( const double  bZ) const
inline

Get the momentum at the start point of the trajectory.

Definition at line 208 of file CDCTrajectory3D.h.

209 {
211 }

◆ getNDF()

size_t getNDF ( ) const
inline

Getter for the number of degrees of freedom of the helix fit.

Definition at line 327 of file CDCTrajectory3D.h.

328 {
329 return getLocalHelix().ndf();
330 }
std::size_t ndf() const
Getter for the number of degrees of freediom used in the helix fit.

◆ getPValue()

double getPValue ( ) const
inline

Getter for p-value.

Definition at line 309 of file CDCTrajectory3D.h.

310 {
311 return TMath::Prob(getChi2(), getNDF());
312 }
double getChi2() const
Getter for the chi2 value of the fit.
size_t getNDF() const
Getter for the number of degrees of freedom of the helix fit.

◆ getSupport()

Vector3D getSupport ( ) const
inline

Getter for the support point of the trajectory in global coordinates, where arcLength2D = 0.

Definition at line 227 of file CDCTrajectory3D.h.

228 {
229 return getLocalHelix()->perigee() + getLocalOrigin();
230 }
Vector3D perigee() const
Getter for the perigee point of the helix.
Definition: Helix.h:234

◆ getTanLambda()

double getTanLambda ( ) const
inline

Getter for the slope of z over the transverse travel distance s.

Definition at line 297 of file CDCTrajectory3D.h.

298 {
299 return getLocalHelix()->tanLambda();
300 }

◆ getTrajectory2D()

CDCTrajectory2D getTrajectory2D ( ) const

Getter for the two dimensional trajectory.

Definition at line 333 of file CDCTrajectory3D.cc.

334{
335 return CDCTrajectory2D(getLocalOrigin().xy(),
336 getLocalHelix().uncertainCircleXY(),
337 getFlightTime());
338}
Particle trajectory as it is seen in xy projection represented as a circle.
double getFlightTime() const
Getter for the time when the particle reached the support point position.

◆ getTrajectorySZ()

CDCTrajectorySZ getTrajectorySZ ( ) const

Getter for the sz trajectory.

Definition at line 340 of file CDCTrajectory3D.cc.

341{
343 globalSZLine.passiveMoveBy(Vector2D(0, -getLocalOrigin().z()));
344 return CDCTrajectorySZ(globalSZLine);
345}
Linear trajectory in sz space.
A line in sz where s is the transverse travel distance as seen in the xy projection with uncertaintie...
void passiveMoveBy(const Vector2D &bySZ)
Moves the coordinate system by the vector by and calculates the new sz line and its covariance matrix...
A two dimensional vector which is equipped with functions for correct handeling of orientation relat...
Definition: Vector2D.h:32

◆ isCurler()

bool isCurler ( double  factor = 1) const

Checks if the trajectory leaves the outer radius of the CDC times the given tolerance factor.

Definition at line 296 of file CDCTrajectory3D.cc.

297{
299 return getMaximalCylindricalR() < factor * topology.getOuterCylindricalR();
300}
double getMaximalCylindricalR() const
Getter for the maximal distance from the origin.
Class representating the sense wire arrangement in the whole of the central drift chamber.
static CDCWireTopology & getInstance()
Getter for the singleton instance of the wire topology.
double getOuterCylindricalR() const
Getter for the outer radius of the outer most wire layer.

◆ isFitted()

bool isFitted ( ) const
inline

Checks if the trajectory has already been set to a valid value.

Definition at line 124 of file CDCTrajectory3D.h.

125 {
126 return not isInvalid();
127 }
bool isInvalid() const
Checks if the trajectory is already set to a valid value.

◆ isInvalid()

bool isInvalid ( ) const
inline

Checks if the trajectory is already set to a valid value.

Definition at line 118 of file CDCTrajectory3D.h.

119 {
120 return m_localHelix->isInvalid();
121 }
bool isInvalid() const
Indicates if the stored parameter combination designates a valid helix.
Definition: Helix.h:74

◆ reconstruct3D()

Vector3D reconstruct3D ( const WireLine wireLine,
double  distance = 0.0 
) const

Gives the three dimensional point which is on the dirft circle away from the wire line.

This method makes the reconstruction of the z coordinate possible by using the skewness
stereo layer of the stereo wires. The point is determined such that it is at the (signed) distance to the wire line.

◆ reverse()

void reverse ( )
inline

Reverses the trajectory in place.

Definition at line 149 of file CDCTrajectory3D.h.

150 {
153 }
void reverse()
Flips the orientation of the circle in place.

◆ reversed()

CDCTrajectory3D reversed ( ) const
inline

Returns the reverse trajectory as a copy.

Definition at line 156 of file CDCTrajectory3D.h.

157 {
159 }
CDCTrajectory3D reversed() const
Returns the reverse trajectory as a copy.

◆ setChi2()

void setChi2 ( const double  chi2)
inline

Setter for the chi square value of the helix fit.

Definition at line 321 of file CDCTrajectory3D.h.

322 {
323 return m_localHelix.setChi2(chi2);
324 }
void setChi2(const double chi2)
Setter for the chi square value of the helix fit.

◆ setFlightTime()

void setFlightTime ( double  flightTime)
inline

Setter for the time when the particle reached the support point position.

Definition at line 370 of file CDCTrajectory3D.h.

371 {
372 m_flightTime = flightTime;
373 }

◆ setLocalHelix()

void setLocalHelix ( const UncertainHelix localHelix)
inline

Setter for the helix that describes the trajectory in local coordinates.

Definition at line 345 of file CDCTrajectory3D.h.

346 {
347 m_localHelix = localHelix;
348 }

◆ setLocalOrigin()

double setLocalOrigin ( const Vector3D localOrigin)

Setter for the origin of the local coordinate system.

This sets the origin point the local helix representation is subjected. The local helix is changed such that the set of points in global space is not changed.

Definition at line 366 of file CDCTrajectory3D.cc.

367{
368 double arcLength2D = calcArcLength2D(localOrigin);
369 double factor2DTo3D = hypot2(1, getTanLambda());
370 double arcLength3D = arcLength2D * factor2DTo3D;
371 m_flightTime += arcLength3D / Const::speedOfLight;
373 m_localOrigin = localOrigin;
374 return arcLength2D;
375}
static const double speedOfLight
[cm/ns]
Definition: Const.h:695
double calcArcLength2D(const Vector3D &point) const
Calculate the travel distance from the start position of the trajectory.
double getTanLambda() const
Getter for the slope of z over the transverse travel distance s.
void passiveMoveBy(const Vector3D &by)
Moves the coordinate system by the vector by and calculates the new perigee and its covariance matrix...

◆ setNDF()

void setNDF ( std::size_t  ndf)
inline

Setter for the number of degrees of freedom of the helix fit.

Definition at line 333 of file CDCTrajectory3D.h.

334 {
335 return m_localHelix.setNDF(ndf);
336 }
void setNDF(std::size_t ndf)
Setter for the number of degrees of freediom used in the helix fit.

◆ shiftPeriod()

double shiftPeriod ( int  nPeriods)

Adjusts the z0 to the one that lies n periods forward.

Returns
The two dimensional arc length needed to travel from the old to the new support point.

Definition at line 326 of file CDCTrajectory3D.cc.

327{
328 double arcLength2D = m_localHelix.shiftPeriod(nPeriods);
329 m_flightTime += arcLength2D / Const::speedOfLight;
330 return arcLength2D;
331}
double shiftPeriod(int nPeriods)
Adjust the arclength measure to start n periods later.

Member Data Documentation

◆ m_flightTime

double m_flightTime = NAN
private

Memory for the estimation of the time at which the particle arrived at the support point.

Definition at line 384 of file CDCTrajectory3D.h.

◆ m_localHelix

UncertainHelix m_localHelix
private

Memory for the generalized circle describing the trajectory in coordinates from the local origin.

Definition at line 381 of file CDCTrajectory3D.h.

◆ m_localOrigin

Vector3D m_localOrigin
private

Memory for local coordinate origin of the circle representing the trajectory in global coordinates.

Definition at line 377 of file CDCTrajectory3D.h.


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