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 CDC::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 46 of file CDCTrajectory3D.h.

Constructor & Destructor Documentation

◆ CDCTrajectory3D() [1/12]

CDCTrajectory3D ( )
inline

Default constructor for ROOT compatibility.

Definition at line 50 of file CDCTrajectory3D.h.

51 : m_localOrigin()
52 , m_localHelix()
53 {
54 }

◆ CDCTrajectory3D() [2/12]

CDCTrajectory3D ( const UncertainHelix & helix)
inlineexplicit

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

Definition at line 57 of file CDCTrajectory3D.h.

58 : m_localOrigin(0.0, 0.0, 0.0)
59 , m_localHelix(helix)
60 , m_flightTime(0.0)
61 {
62 }

◆ CDCTrajectory3D() [3/12]

CDCTrajectory3D ( const Helix & helix)
inlineexplicit

conversion constructor to make that one stupid test work

Definition at line 65 of file CDCTrajectory3D.h.

65 : CDCTrajectory3D(UncertainHelix(helix))
66 {
67 }

◆ 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 70 of file CDCTrajectory3D.h.

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

◆ 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 48 of file CDCTrajectory3D.cc.

50 : m_localOrigin(trajectory2D.getLocalOrigin())
51 , m_localHelix(trajectory2D.getLocalCircle(), trajectorySZ.getSZLine())
52 , m_flightTime(trajectory2D.getFlightTime())
53{
54}
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 circle in local coordinates.
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
Memory for the estimation of the time at which the particle arrived at the support point.
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 56 of file CDCTrajectory3D.cc.

58{
59}
CDCTrajectory3D()
Default constructor for ROOT compatibility.
static CDCTrajectorySZ basicAssumption()
Constructs a basic assumption, what the z0 start position and the sz slope are, including some broad ...

◆ 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 61 of file CDCTrajectory3D.cc.

66 : m_localOrigin(pos3D)
68 mom3D.xy().unit(),
69 0.0,
70 mom3D.cotTheta(),
71 0.0)
72 , m_flightTime(time)
73{
74}
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:339
double norm() const
Calculates the length of the vector.
Definition Vector2D.h:193
const Vector2D & xy() const
Getter for the xy projected vector ( reference ! )
Definition Vector3D.h:511
double cotTheta() const
Getter for the cotangent of the polar angle.
Definition Vector3D.h:561

◆ 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 76 of file CDCTrajectory3D.cc.

80 : CDCTrajectory3D(pos3D, time, mom3D, charge, CDCBFieldUtil::getBFieldZ(pos3D))
81{
82}
static double getBFieldZ()
Getter for the signed magnetic field strength 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 84 of file CDCTrajectory3D.cc.

86 mcParticle.getProductionTime(),
87 Vector3D{mcParticle.getMomentum()},
88 mcParticle.getCharge(),
89 bZ)
90{
91}
ROOT::Math::XYZVector getProductionVertex() const
Return production vertex position.
Definition MCParticle.h:178
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:148
ROOT::Math::XYZVector getMomentum() const
Return momentum.
Definition MCParticle.h:187
HepGeom::Vector3D< double > Vector3D
3D Vector
Definition Cell.h:34

◆ CDCTrajectory3D() [10/12]

CDCTrajectory3D ( const MCParticle & mcParticle)
explicit

Construct a trajectory from the MCParticles vertex and momentum.

Definition at line 93 of file CDCTrajectory3D.cc.

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

◆ 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 101 of file CDCTrajectory3D.cc.

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

◆ CDCTrajectory3D() [12/12]

CDCTrajectory3D ( const genfit::TrackCand & gfTrackCand)
explicit

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

Definition at line 171 of file CDCTrajectory3D.cc.

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

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 179 of file CDCTrajectory3D.h.

180 {
181 return getLocalHelix()->circleXY().arcLengthTo((point - getLocalOrigin()).xy());
182 }

◆ clear()

void clear ( )
inline

Clears all information from this trajectoy.

Definition at line 131 of file CDCTrajectory3D.h.

132 {
133 m_localOrigin.set(0.0, 0.0, 0.0);
134 m_localHelix.invalidate();
135 m_flightTime = NAN;
136 }

◆ fillInto() [1/2]

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

Copies the trajectory information to the Genfit track candidate.

Definition at line 267 of file CDCTrajectory3D.cc.

268{
269 // Set the start parameters
270 Vector3D position = getSupport();
271 Vector3D momentum = bZ == 0 ? getFlightDirection3DAtSupport() : getMom3DAtSupport(bZ);
272 ESign charge = getChargeSign();
273
274 // Do not propagate invalid fits, signal that the fit is invalid to the caller.
275 if (not ESignUtil::isValid(charge) or momentum.hasNAN() or position.hasNAN()) {
276 return false;
277 }
278
279 gfTrackCand.setPosMomSeed(position, momentum, charge);
280
281 const CovarianceMatrix<6> cov6 = calculateCovarianceMatrix(getLocalHelix(), momentum, charge, bZ);
282 TMatrixDSym covSeed = TMatrixConversion::toTMatrix(cov6);
283
284 gfTrackCand.setCovSeed(covSeed);
285
286 return true;
287}
Vector3D getMom3DAtSupport(const double bZ) const
Get the momentum at the start point of the trajectory.
ESign getChargeSign() const
Gets the charge sign of the trajectory.
Vector3D getFlightDirection3DAtSupport() const
Get the unit momentum at the start point of the trajectory.
const UncertainHelix & getLocalHelix() const
Getter for the helix in local coordinates.
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 261 of file CDCTrajectory3D.cc.

262{
263 Vector3D position = getSupport();
264 return fillInto(trackCand, CDCBFieldUtil::getBFieldZ(position));
265}
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 317 of file CDCTrajectory3D.cc.

318{
319 Vector3D position = getSupport();
320 double tanLambda = getLocalHelix()->tanLambda();
321 double factor2DTo3D = hypot2(1, tanLambda);
322 double curvatureXY = getLocalHelix()->curvatureXY();
323 double absMom2D = CDCBFieldUtil::curvatureToAbsMom2D(curvatureXY, position);
324 return factor2DTo3D * absMom2D;
325}
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 308 of file CDCTrajectory3D.cc.

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

◆ getArcLength2DPeriod()

double getArcLength2DPeriod ( ) const
inline

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

Definition at line 185 of file CDCTrajectory3D.h.

186 {
187 return getLocalHelix()->arcLength2DPeriod();
188 }

◆ getCartesianCovariance()

CovarianceMatrix< 6 > getCartesianCovariance ( double bZ) const

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

Definition at line 289 of file CDCTrajectory3D.cc.

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

◆ getChargeSign()

ESign getChargeSign ( ) const

Gets the charge sign of the trajectory.

Definition at line 303 of file CDCTrajectory3D.cc.

304{
305 return CDCBFieldUtil::ccwInfoToChargeSign(getLocalHelix()->circleXY().orientation());
306}
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 316 of file CDCTrajectory3D.h.

317 {
318 return getLocalHelix().chi2();
319 }

◆ getCurvatureXY()

double getCurvatureXY ( ) const
inline

Getter for the curvature as seen from the xy projection.

Definition at line 304 of file CDCTrajectory3D.h.

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

◆ getFlightDirection3DAtSupport()

Vector3D getFlightDirection3DAtSupport ( ) const
inline

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

Definition at line 221 of file CDCTrajectory3D.h.

222 {
223 return getLocalHelix()->tangential();
224 }

◆ getFlightTime()

double getFlightTime ( ) const
inline

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

Definition at line 365 of file CDCTrajectory3D.h.

366 {
367 return m_flightTime;
368 }

◆ getGlobalCenter()

Vector2D getGlobalCenter ( ) const
inline

Getter for the center of the helix in global coordinates.

Definition at line 240 of file CDCTrajectory3D.h.

241 {
242 return getLocalHelix()->centerXY() + m_localOrigin.xy();
243 }

◆ getGlobalCircle()

PerigeeCircle getGlobalCircle ( ) const

Getter for the circle in global coordinates.

Definition at line 349 of file CDCTrajectory3D.cc.

350{
351 // Down cast since we do not necessarily won't the covariance matrix transformed as well
352 PerigeeCircle result(getLocalHelix()->circleXY());
353 result.passiveMoveBy(-getLocalOrigin().xy());
354 return result;
355}
const Vector3D & getLocalOrigin() const
Getter for the origin of the local coordinate system.

◆ getGlobalImpact()

double getGlobalImpact ( ) const
inline

Getter for the signed impact parameter of the trajectory.

Definition at line 263 of file CDCTrajectory3D.h.

264 {
265 return getLocalHelix()->distanceXY(-m_localOrigin.xy());
266 }

◆ getGlobalPerigee()

Vector3D getGlobalPerigee ( ) const
inline

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

Definition at line 234 of file CDCTrajectory3D.h.

235 {
236 return getLocalHelix()->closestXY(-m_localOrigin.xy()) + m_localOrigin;
237 }

◆ getLocalCircle()

UncertainPerigeeCircle getLocalCircle ( ) const

Getter for the circle in local coordinates.

Definition at line 357 of file CDCTrajectory3D.cc.

358{
360}
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 285 of file CDCTrajectory3D.h.

286 {
287 return getLocalHelix().covariance(iRow, iCol);
288 }

◆ getLocalHelix()

const UncertainHelix & getLocalHelix ( ) const
inline

Getter for the helix in local coordinates.

Definition at line 340 of file CDCTrajectory3D.h.

341 {
342 return m_localHelix;
343 }

◆ getLocalOrigin()

const Vector3D & getLocalOrigin ( ) const
inline

Getter for the origin of the local coordinate system.

Definition at line 352 of file CDCTrajectory3D.h.

353 {
354 return m_localOrigin;
355 }

◆ getLocalSZLine()

UncertainSZLine getLocalSZLine ( ) const

Getter for the sz line starting from the local origin.

Definition at line 362 of file CDCTrajectory3D.cc.

363{
365}
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 292 of file CDCTrajectory3D.h.

293 {
294 return getLocalHelix().variance(i);
295 }

◆ getMaximalCylindricalR()

double getMaximalCylindricalR ( ) const
inline

Getter for the maximal distance from the origin.

Definition at line 251 of file CDCTrajectory3D.h.

252 {
253 return std::fabs(getGlobalImpact() + 2 * getLocalHelix()->radiusXY());
254 }

◆ getMinimalCylindricalR()

double getMinimalCylindricalR ( ) const
inline

Getter for the minimal distance from the origin.

Definition at line 257 of file CDCTrajectory3D.h.

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

◆ getMom3DAtSupport() [1/2]

Vector3D getMom3DAtSupport ( ) const
inline

Get the momentum at the start point of the trajectory.

Definition at line 215 of file CDCTrajectory3D.h.

216 {
217 return getFlightDirection3DAtSupport() *= getAbsMom3D();
218 }

◆ getMom3DAtSupport() [2/2]

Vector3D getMom3DAtSupport ( const double bZ) const
inline

Get the momentum at the start point of the trajectory.

Definition at line 209 of file CDCTrajectory3D.h.

210 {
211 return getFlightDirection3DAtSupport() *= getAbsMom3D(bZ);
212 }

◆ getNDF()

size_t getNDF ( ) const
inline

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

Definition at line 328 of file CDCTrajectory3D.h.

329 {
330 return getLocalHelix().ndf();
331 }

◆ getPValue()

double getPValue ( ) const
inline

Getter for p-value.

Definition at line 310 of file CDCTrajectory3D.h.

311 {
312 return TMath::Prob(getChi2(), getNDF());
313 }

◆ getSupport()

Vector3D getSupport ( ) const
inline

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

Definition at line 228 of file CDCTrajectory3D.h.

229 {
230 return getLocalHelix()->perigee() + getLocalOrigin();
231 }

◆ getTanLambda()

double getTanLambda ( ) const
inline

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

Definition at line 298 of file CDCTrajectory3D.h.

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

◆ getTrajectory2D()

CDCTrajectory2D getTrajectory2D ( ) const

Getter for the two dimensional trajectory.

Definition at line 334 of file CDCTrajectory3D.cc.

335{
336 return CDCTrajectory2D(getLocalOrigin().xy(),
337 getLocalHelix().uncertainCircleXY(),
338 getFlightTime());
339}
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 341 of file CDCTrajectory3D.cc.

342{
343 UncertainSZLine globalSZLine = getLocalHelix().uncertainSZLine();
344 globalSZLine.passiveMoveBy(Vector2D(0, -getLocalOrigin().z()));
345 return CDCTrajectorySZ(globalSZLine);
346}
void passiveMoveBy(const Vector2D &bySZ)
Moves the coordinate system by the vector by and calculates the new sz line and its covariance matrix...

◆ 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 297 of file CDCTrajectory3D.cc.

298{
299 const CDCWireTopology& topology = CDCWireTopology::getInstance();
300 return getMaximalCylindricalR() < factor * topology.getOuterCylindricalR();
301}
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.
double getMaximalCylindricalR() const
Getter for the maximal distance from the origin.

◆ isFitted()

bool isFitted ( ) const
inline

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

Definition at line 125 of file CDCTrajectory3D.h.

126 {
127 return not isInvalid();
128 }

◆ isInvalid()

bool isInvalid ( ) const
inline

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

Definition at line 119 of file CDCTrajectory3D.h.

120 {
121 return m_localHelix->isInvalid();
122 }

◆ reconstruct3D()

Vector3D reconstruct3D ( const CDC::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 150 of file CDCTrajectory3D.h.

151 {
152 m_localHelix.reverse();
153 m_flightTime = -m_flightTime;
154 }

◆ reversed()

CDCTrajectory3D reversed ( ) const
inline

Returns the reverse trajectory as a copy.

Definition at line 157 of file CDCTrajectory3D.h.

158 {
159 return CDCTrajectory3D(getLocalOrigin(), getLocalHelix().reversed(), -getFlightTime());
160 }
EForwardBackward reversed(EForwardBackward eForwardBackward)
Return the reversed forward backward indicator. Leaves EForwardBackward::c_Invalid the same.

◆ setChi2()

void setChi2 ( const double chi2)
inline

Setter for the chi square value of the helix fit.

Definition at line 322 of file CDCTrajectory3D.h.

323 {
324 return m_localHelix.setChi2(chi2);
325 }

◆ setFlightTime()

void setFlightTime ( double flightTime)
inline

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

Definition at line 371 of file CDCTrajectory3D.h.

372 {
373 m_flightTime = flightTime;
374 }

◆ setLocalHelix()

void setLocalHelix ( const UncertainHelix & localHelix)
inline

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

Definition at line 346 of file CDCTrajectory3D.h.

347 {
348 m_localHelix = localHelix;
349 }

◆ 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 367 of file CDCTrajectory3D.cc.

368{
369 double arcLength2D = calcArcLength2D(localOrigin);
370 double factor2DTo3D = hypot2(1, getTanLambda());
371 double arcLength3D = arcLength2D * factor2DTo3D;
372 m_flightTime += arcLength3D / Const::speedOfLight;
373 m_localHelix.passiveMoveBy(localOrigin - m_localOrigin);
374 m_localOrigin = localOrigin;
375 return arcLength2D;
376}
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.

◆ setNDF()

void setNDF ( std::size_t ndf)
inline

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

Definition at line 334 of file CDCTrajectory3D.h.

335 {
336 return m_localHelix.setNDF(ndf);
337 }

◆ 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 327 of file CDCTrajectory3D.cc.

328{
329 double arcLength2D = m_localHelix.shiftPeriod(nPeriods);
330 m_flightTime += arcLength2D / Const::speedOfLight;
331 return arcLength2D;
332}

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 385 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 382 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 378 of file CDCTrajectory3D.h.


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