Belle II Software development

This class represents an ideal helix in perigee parameterization. More...

#include <Helix.h>

Inheritance diagram for Helix:
RelationsInterface< BASE >

Public Member Functions

 Helix ()
 Constructor initializing all perigee parameters to zero.
 
 Helix (const ROOT::Math::XYZVector &position, const ROOT::Math::XYZVector &momentum, const short int charge, const double bZ)
 Constructor initializing class with a fit result.
 
 Helix (const double &d0, const double &phi0, const double &omega, const double &z0, const double &tanLambda)
 Constructor initializing class with perigee parameters.
 
void addRelationTo (const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
 Add a relation from this object to another object (with caching).
 
void addRelationTo (const TObject *object, float weight=1.0, const std::string &namedRelation="") const
 Add a relation from this object to another object (no caching, can be quite slow).
 
void copyRelations (const RelationsInterface< BASE > *sourceObj)
 Copies all relations of sourceObj (pointing from or to sourceObj) to this object (including weights).
 
template<class TO >
RelationVector< TO > getRelationsTo (const std::string &name="", const std::string &namedRelation="") const
 Get the relations that point from this object to another store array.
 
template<class FROM >
RelationVector< FROM > getRelationsFrom (const std::string &name="", const std::string &namedRelation="") const
 Get the relations that point from another store array to this object.
 
template<class T >
RelationVector< T > getRelationsWith (const std::string &name="", const std::string &namedRelation="") const
 Get the relations between this object and another store array.
 
template<class TO >
TO * getRelatedTo (const std::string &name="", const std::string &namedRelation="") const
 Get the object to which this object has a relation.
 
template<class FROM >
FROM * getRelatedFrom (const std::string &name="", const std::string &namedRelation="") const
 Get the object from which this object has a relation.
 
template<class T >
T * getRelated (const std::string &name="", const std::string &namedRelation="") const
 Get the object to or from which this object has a relation.
 
template<class TO >
std::pair< TO *, float > getRelatedToWithWeight (const std::string &name="", const std::string &namedRelation="") const
 Get first related object & weight of relation pointing to an array.
 
template<class FROM >
std::pair< FROM *, float > getRelatedFromWithWeight (const std::string &name="", const std::string &namedRelation="") const
 Get first related object & weight of relation pointing from an array.
 
template<class T >
std::pair< T *, float > getRelatedWithWeight (const std::string &name="", const std::string &namedRelation="") const
 Get first related object & weight of relation pointing from/to an array.
 
virtual std::string getName () const
 Return a short name that describes this object, e.g.
 
virtual std::string getInfoHTML () const
 Return a short summary of this object's contents in HTML format.
 
std::string getInfo () const
 Return a short summary of this object's contents in raw text format.
 
std::string getArrayName () const
 Get name of array this object is stored in, or "" if not found.
 
int getArrayIndex () const
 Returns this object's array index (in StoreArray), or -1 if not found.
 
Getters for perigee helix parameters
double getD0 () const
 Getter for d0, which is the signed distance to the perigee in the r-phi plane.
 
double getPhi0 () const
 Getter for phi0, which is the azimuth angle of the transverse momentum at the perigee.
 
double getCosPhi0 () const
 Getter for the cosine of the azimuth angle of travel direction at the perigee.
 
double getSinPhi0 () const
 Getter for the cosine of the azimuth angle of travel direction at the perigee.
 
double getOmega () const
 Getter for omega, which is a signed curvature measure of the track.
 
double getZ0 () const
 Getter for z0, which is the z coordinate of the perigee.
 
double getTanLambda () const
 Getter for tan lambda, which is the z over two dimensional arc length slope of the track.
 
double getCotTheta () const
 Getter for cot theta, which is the z over two dimensional arc length slope of the track.
 

Protected Member Functions

void calcArcLength2DAndDrAtXY (const double &x, const double &y, double &arcLength2D, double &dr) const
 Helper method to calculate the signed two dimensional arc length and the signed distance to the circle of a point in the xy projection.
 
double calcArcLength2DAtDeltaCylindricalRAndDr (const double &deltaCylindricalR, const double &dr) const
 Helper method to calculate the two dimensional arc length from the perigee to a point at cylindrical radius, which also has the distance dr from the circle in the xy projection.
 
TClonesArray * getArrayPointer () const
 Returns the pointer to the raw DataStore array holding this object (protected since these arrays are easy to misuse).
 

Static Protected Member Functions

static double calcASinXDividedByX (const double &x)
 Implementation of the function asin(x) / x which handles small x values smoothly.
 
static double calcATanXDividedByX (const double &x)
 Implementation of the function atan(x) / x which handles small x values smoothly.
 
static double calcDerivativeOfATanXDividedByX (const double &x)
 Implementation of the function d / dx (atan(x) / x) which handles small x values smoothly.
 

Private Member Functions

void setCartesian (const ROOT::Math::XYZVector &position, const ROOT::Math::XYZVector &momentum, const short int charge, const double bZ)
 Cartesian to Perigee conversion.
 
 ClassDef (Helix, 2)
 This class represents an ideal helix in perigee parameterization.
 

Private Attributes

Double32_t m_d0
 Memory for the signed distance to the perigee.
 
Double32_t m_phi0
 Memory for the azimuth angle between the transverse momentum and the x axis, which is in [-pi, pi].
 
Double32_t m_omega
 Memory for the curvature of the signed curvature.
 
Double32_t m_z0
 Memory for the z coordinate of the perigee.
 
Double32_t m_tanLambda
 Memory for the slope of the track in the z coordinate over the two dimensional arc length (dz/ds)
 
DataStore::StoreEntrym_cacheDataStoreEntry
 Cache of the data store entry to which this object belongs.
 
int m_cacheArrayIndex
 Cache of the index in the TClonesArray to which this object belongs.
 

Friends

std::ostream & operator<< (std::ostream &output, const Helix &helix)
 Output operator for debugging and the generation of unittest error messages.
 

Getters for cartesian parameters of the perigee

double getPerigeeX () const
 Calculates the x coordinate of the perigee point.
 
double getPerigeeY () const
 Calculates the y coordinate of the perigee point.
 
double getPerigeeZ () const
 Calculates the z coordinate of the perigee point.
 
ROOT::Math::XYZVector getPerigee () const
 Getter for the perigee position.
 
double getMomentumX (const double bZ) const
 Calculates the x momentum of the particle at the perigee point.
 
double getMomentumY (const double bZ) const
 Calculates the y momentum of the particle at the perigee point.
 
double getMomentumZ (const double bZ) const
 Calculates the z momentum of the particle at the perigee point.
 
ROOT::Math::XYZVector getMomentum (const double bZ) const
 Getter for vector of momentum at the perigee position.
 
ROOT::Math::XYZVector getDirection () const
 Getter for unit vector of momentum at the perigee position.
 
double getTransverseMomentum (const double bZ) const
 Getter for the absolute value of the transverse momentum at the perigee.
 
double getKappa (const double bZ) const
 Getter for kappa, which is charge / transverse momentum or equivalently omega * alpha.
 
short getChargeSign () const
 Return track charge sign (1, 0 or -1).
 
static double getAlpha (const double bZ)
 Calculates the alpha value for a given magnetic field in Tesla.
 

Simple extrapolations of the ideal helix

double getArcLength2DAtCylindricalR (const double &cylindricalR) const
 Calculates the transverse travel distance at the point the helix first reaches the given cylindrical radius.
 
double getArcLength2DAtXY (const double &x, const double &y) const
 Calculates the two dimensional arc length at which the circle in the xy projection is closest to the point.
 
double getArcLength2DAtNormalPlane (const double &x, const double &y, const double &nX, const double &nY) const
 Calculates the arc length to reach a plane parallel to the z axes.
 
ROOT::Math::XYZVector getPositionAtArcLength2D (const double &arcLength2D) const
 Calculates the position on the helix at the given two dimensional arc length.
 
ROOT::Math::XYZVector getTangentialAtArcLength2D (const double &arcLength2D) const
 Calculates the tangential vector to the helix curve at the given two dimensional arc length.
 
ROOT::Math::XYZVector getUnitTangentialAtArcLength2D (const double &arcLength2D) const
 Calculates the unit tangential vector to the helix curve at the given two dimensional arc length.
 
ROOT::Math::XYZVector getMomentumAtArcLength2D (const double &arcLength2D, const double &bz) const
 Calculates the momentum vector at the given two dimensional arc length.
 
double passiveMoveBy (const ROOT::Math::XYZVector &by)
 Moves origin of the coordinate system (passive transformation) by the given vector.
 
double passiveMoveBy (const double &byX, const double &byY, const double &byZ)
 Moves origin of the coordinate system (passive transformation) orthogonal to the z axis by the given vector.
 
TMatrixD calcPassiveMoveByJacobian (const ROOT::Math::XYZVector &by, const double expandBelowChi=M_PI/8) const
 Calculate the 5x5 jacobian matrix for the transport of the helix parameters, when moving the origin of the coordinate system to a new location.
 
void calcPassiveMoveByJacobian (const double &byX, const double &byY, TMatrixD &jacobian, const double expandBelowChi=M_PI/8) const
 Calculate the jacobian matrix for the transport of the helix parameters, when moving the origin of the coordinate system to a new location.
 
void reverse ()
 Reverses the direction of travel of the helix in place.
 
double calcArcLength2DFromSecantLength (const double &secantLength2D) const
 Helper function to calculate the two dimensional arc length from the length of a secant.
 
double calcSecantLengthToArcLength2DFactor (const double &secantLength2D) const
 Helper function to calculate the factor between the dimensional secant length and the two dimensional arc length as seen in xy projection of the helix.
 
static double reversePhi (const double &phi)
 Reverses an azimuthal angle to the opposite direction.
 

Detailed Description

This class represents an ideal helix in perigee parameterization.

The used perigee parameters are:

  1. $ d_0 $ - the signed distance from the origin to the perigee. The sign is positive (negative), if the angle from the xy perigee position vector to the transverse momentum vector is +pi/2 (-pi/2). $d_0$ has the same sign as getPerigee().Cross(getMomentum()).Z().
  2. $ \phi_0 $ - the angle in the xy projection between the transverse momentum and the x axis, which is in [-pi, pi]
  3. $ \omega $ - the signed curvature of the track where the sign is given by the charge of the particle
  4. $ z_0 $ - z coordinate of the perigee
  5. $ \tan \lambda $ - the slope of the track in the sz plane (dz/ds)

in that exact order.

Each point on the helix can be adressed by the two dimensional arc length s, which has to be traversed to get to it from the perigee. More precisely the two dimensional arc length means the transverse part of the particles travel distance, hence the arc length of the circle in the xy projection.

If you need different kind of methods / interfaces to the helix please do not hesitate to contact olive.nosp@m.r.fr.nosp@m.ost@d.nosp@m.esy..nosp@m.de Contributions are always welcome.

Definition at line 59 of file Helix.h.

Constructor & Destructor Documentation

◆ Helix() [1/2]

Helix ( const ROOT::Math::XYZVector &  position,
const ROOT::Math::XYZVector &  momentum,
const short int  charge,
const double  bZ 
)

Constructor initializing class with a fit result.

The given position and momentum are extrapolated to the perigee assuming a homogeneous magnetic field in the z direction.

Parameters
positionPosition of the track at the perigee.
momentumMomentum of the track at the perigee.
chargeCharge of the particle.
bZMagnetic field to be used for the calculation of the curvature. It is assumed, that the B-field is homogeneous parallel to the z axis.

Definition at line 33 of file Helix.cc.

37{
38 setCartesian(position, momentum, charge, bZ);
39}
Hep3Vector momentum(double dPhi=0.) const
returns momentum vector after rotating angle dPhi in phi direction.
Definition: Helix.cc:259
void setCartesian(const ROOT::Math::XYZVector &position, const ROOT::Math::XYZVector &momentum, const short int charge, const double bZ)
Cartesian to Perigee conversion.
Definition: Helix.cc:588

◆ Helix() [2/2]

Helix ( const double &  d0,
const double &  phi0,
const double &  omega,
const double &  z0,
const double &  tanLambda 
)

Constructor initializing class with perigee parameters.

Parameters
d0The signed distance from the origin to the perigee. The sign is positive (negative), if the angle from the xy perigee position vector to the transverse momentum vector is +pi/2 (-pi/2). d0 has the same sign as getPosition().Cross(getMomentum()).Z().
phi0The angle between the transverse momentum and the x axis and in [-pi, pi].
omegaThe signed curvature of the track where the sign is given by the charge of the particle.
z0The z coordinate of the perigee.
tanLambdaThe slope of the track in the sz plane (dz/ds).

Definition at line 41 of file Helix.cc.

46 : m_d0(d0),
47 m_phi0(phi0),
48 m_omega(omega),
49 m_z0(z0),
50 m_tanLambda(tanLambda)
51{
52}
double phi0(void) const
Return helix parameter phi0.
Definition: Helix.h:388
Double32_t m_tanLambda
Memory for the slope of the track in the z coordinate over the two dimensional arc length (dz/ds)
Definition: Helix.h:426
Double32_t m_omega
Memory for the curvature of the signed curvature.
Definition: Helix.h:420
Double32_t m_z0
Memory for the z coordinate of the perigee.
Definition: Helix.h:423
Double32_t m_d0
Memory for the signed distance to the perigee.
Definition: Helix.h:414
Double32_t m_phi0
Memory for the azimuth angle between the transverse momentum and the x axis, which is in [-pi,...
Definition: Helix.h:417

Member Function Documentation

◆ addRelationTo() [1/2]

void addRelationTo ( const RelationsInterface< BASE > *  object,
float  weight = 1.0,
const std::string &  namedRelation = "" 
) const
inlineinherited

Add a relation from this object to another object (with caching).

Parameters
objectThe object to which the relation should point.
weightThe weight of the relation.
namedRelationAdditional name for the relation, or "" for the default naming

Definition at line 142 of file RelationsObject.h.

143 {
144 if (object)
146 object, object->m_cacheDataStoreEntry, object->m_cacheArrayIndex, weight, namedRelation);
147 }
void addRelation(const TObject *fromObject, StoreEntry *&fromEntry, int &fromIndex, const TObject *toObject, StoreEntry *&toEntry, int &toIndex, float weight, const std::string &namedRelation)
Add a relation from an object in a store array to another object in a store array.
Definition: DataStore.cc:492
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:54
DataStore::StoreEntry * m_cacheDataStoreEntry
Cache of the data store entry to which this object belongs.
int m_cacheArrayIndex
Cache of the index in the TClonesArray to which this object belongs.

◆ addRelationTo() [2/2]

void addRelationTo ( const TObject *  object,
float  weight = 1.0,
const std::string &  namedRelation = "" 
) const
inlineinherited

Add a relation from this object to another object (no caching, can be quite slow).

Parameters
objectThe object to which the relation should point.
weightThe weight of the relation.
namedRelationAdditional name for the relation, or "" for the default naming

Definition at line 155 of file RelationsObject.h.

156 {
157 StoreEntry* toEntry = nullptr;
158 int toIndex = -1;
159 DataStore::Instance().addRelation(this, m_cacheDataStoreEntry, m_cacheArrayIndex, object, toEntry, toIndex, weight, namedRelation);
160 }

◆ calcArcLength2DAndDrAtXY()

void calcArcLength2DAndDrAtXY ( const double &  x,
const double &  y,
double &  arcLength2D,
double &  dr 
) const
protected

Helper method to calculate the signed two dimensional arc length and the signed distance to the circle of a point in the xy projection.

This function is an implementation detail that prevents some code duplication.

Parameters
xX coordinate of the point to which to extrapolate
yY coordinate of the point to which to extrapolate
[out]arcLength2DThe two dimensional arc length from the perigee at which the closest approach is reached
[out]drSigned distance of the point to circle in the xy projection.

Definition at line 545 of file Helix.cc.

546{
547 // Prepare common variables
548 const double omega = getOmega();
549 const double cosPhi0 = getCosPhi0();
550 const double sinPhi0 = getSinPhi0();
551 const double d0 = getD0();
552
553 const double deltaParallel = x * cosPhi0 + y * sinPhi0;
554 const double deltaOrthogonal = y * cosPhi0 - x * sinPhi0 + d0;
555 const double deltaCylindricalR = hypot(deltaOrthogonal, deltaParallel);
556 const double deltaCylindricalRSquared = deltaCylindricalR * deltaCylindricalR;
557
558 const double A = 2 * deltaOrthogonal + omega * deltaCylindricalRSquared;
559 const double U = sqrt(1 + omega * A);
560 const double UOrthogonal = 1 + omega * deltaOrthogonal; // called nu in the Karimaki paper.
561 const double UParallel = omega * deltaParallel;
562 // Note U is a vector pointing from the middle of the projected circle scaled by a factor omega.
563
564 // Calculate dr
565 dr = A / (1 + U);
566
567 // Calculate the absolute value of the arc length
568 const double chi = -atan2(UParallel, UOrthogonal);
569
570 if (fabs(chi) < M_PI / 8) { // Rough guess where the critical zone for approximations begins
571 // Close side of the circle
572 double principleArcLength2D = deltaParallel / UOrthogonal;
573 arcLength2D = principleArcLength2D * calcATanXDividedByX(principleArcLength2D * omega);
574 } else {
575 // Far side of the circle
576 // If the far side of the circle is a well definied concept meaning that we have big enough omega.
577 arcLength2D = -chi / omega;
578 }
579}
double dr(void) const
Return helix parameter dr.
Definition: Helix.h:380
HepPoint3D x(double dPhi=0.) const
returns position after rotating angle dPhi in phi direction.
Definition: Helix.cc:197
double cosPhi0(void) const
Return cos phi0.
Definition: Helix.h:500
double sinPhi0(void) const
Return sin phi0.
Definition: Helix.h:492
double getSinPhi0() const
Getter for the cosine of the azimuth angle of travel direction at the perigee.
Definition: Helix.h:384
double getCosPhi0() const
Getter for the cosine of the azimuth angle of travel direction at the perigee.
Definition: Helix.h:381
double getOmega() const
Getter for omega, which is a signed curvature measure of the track.
Definition: Helix.h:387
double getD0() const
Getter for d0, which is the signed distance to the perigee in the r-phi plane.
Definition: Helix.h:372
static double calcATanXDividedByX(const double &x)
Implementation of the function atan(x) / x which handles small x values smoothly.
Definition: Helix.cc:477
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28

◆ calcArcLength2DAtDeltaCylindricalRAndDr()

double calcArcLength2DAtDeltaCylindricalRAndDr ( const double &  deltaCylindricalR,
const double &  dr 
) const
protected

Helper method to calculate the two dimensional arc length from the perigee to a point at cylindrical radius, which also has the distance dr from the circle in the xy projection.

This function is an implementation detail that prevents some code duplication.

Parameters
deltaCylindricalRThe absolute distance of the point in question to the perigee in the xy projection
drSigned distance of the point to circle in the xy projection.
Returns
The absolute two dimensional arc length from the perigee to the point.

Definition at line 581 of file Helix.cc.

582{
583 const double omega = getOmega();
584 double secantLength = sqrt((deltaCylindricalR + dr) * (deltaCylindricalR - dr) / (1 + dr * omega));
585 return calcArcLength2DFromSecantLength(secantLength);
586}
double calcArcLength2DFromSecantLength(const double &secantLength2D) const
Helper function to calculate the two dimensional arc length from the length of a secant.
Definition: Helix.cc:426

◆ calcArcLength2DFromSecantLength()

double calcArcLength2DFromSecantLength ( const double &  secantLength2D) const

Helper function to calculate the two dimensional arc length from the length of a secant.

Translates the direct length between two point on the circle in the xy projection to the two dimensional arc length on the circle Behaves smoothly in the limit of vanishing curvature.

Definition at line 426 of file Helix.cc.

427{
428 return secantLength * calcSecantLengthToArcLength2DFactor(secantLength);
429}
double calcSecantLengthToArcLength2DFactor(const double &secantLength2D) const
Helper function to calculate the factor between the dimensional secant length and the two dimensional...
Definition: Helix.cc:432

◆ calcASinXDividedByX()

double calcASinXDividedByX ( const double &  x)
staticprotected

Implementation of the function asin(x) / x which handles small x values smoothly.

Definition at line 438 of file Helix.cc.

439{
440 // Approximation of asin(x) / x
441 // Inspired by BOOST's sinc
442
443 BOOST_MATH_STD_USING;
444
445 auto const taylor_n_bound = boost::math::tools::forth_root_epsilon<double>();
446
447 if (abs(x) >= taylor_n_bound) {
448 if (fabs(x) == 1) {
449 return M_PI / 2.0;
450
451 } else {
452 return asin(x) / x;
453
454 }
455
456 } else {
457 // approximation by taylor series in x at 0 up to order 0
458 double result = 1.0;
459
460 auto const taylor_0_bound = boost::math::tools::epsilon<double>();
461 if (abs(x) >= taylor_0_bound) {
462 double x2 = x * x;
463 // approximation by taylor series in x at 0 up to order 2
464 result += x2 / 6.0;
465
466 auto const taylor_2_bound = boost::math::tools::root_epsilon<double>();
467 if (abs(x) >= taylor_2_bound) {
468 // approximation by taylor series in x at 0 up to order 4
469 result += x2 * x2 * (3.0 / 40.0);
470 }
471 }
472 return result;
473 }
474
475}

◆ calcATanXDividedByX()

double calcATanXDividedByX ( const double &  x)
staticprotected

Implementation of the function atan(x) / x which handles small x values smoothly.

Definition at line 477 of file Helix.cc.

478{
479 // Approximation of atan(x) / x
480 // Inspired by BOOST's sinc
481
482 BOOST_MATH_STD_USING;
483
484 auto const taylor_n_bound = boost::math::tools::forth_root_epsilon<double>();
485
486 if (abs(x) >= taylor_n_bound) {
487 return atan(x) / x;
488
489 } else {
490 // approximation by taylor series in x at 0 up to order 0
491 double result = 1.0;
492
493 auto const taylor_0_bound = boost::math::tools::epsilon<double>();
494 if (abs(x) >= taylor_0_bound) {
495 double x2 = x * x;
496 // approximation by taylor series in x at 0 up to order 2
497 result -= x2 / 3.0;
498
499 auto const taylor_2_bound = boost::math::tools::root_epsilon<double>();
500 if (abs(x) >= taylor_2_bound) {
501 // approximation by taylor series in x at 0 up to order 4
502 result += x2 * x2 * (1.0 / 5.0);
503 }
504 }
505 return result;
506 }
507}
double atan(double a)
atan for double
Definition: beamHelpers.h:34

◆ calcDerivativeOfATanXDividedByX()

double calcDerivativeOfATanXDividedByX ( const double &  x)
staticprotected

Implementation of the function d / dx (atan(x) / x) which handles small x values smoothly.

Definition at line 509 of file Helix.cc.

510{
511 // Approximation of atan(x) / x
512 // Inspired by BOOST's sinc
513
514 BOOST_MATH_STD_USING;
515
516 auto const taylor_n_bound = boost::math::tools::forth_root_epsilon<double>();
517
518 const double x2 = x * x;
519 if (abs(x) >= taylor_n_bound) {
520 const double chi = atan(x);
521 return ((1 - chi / x) / x - chi) / (1 + x2);
522
523 } else {
524 // approximation by taylor series in x at 0 up to order 0
525 double result = 1.0;
526
527 auto const taylor_0_bound = boost::math::tools::epsilon<double>();
528 if (abs(x) >= taylor_0_bound) {
529 // approximation by taylor series in x at 0 up to order 2
530 result -= 2.0 * x / 3.0;
531
532 auto const taylor_2_bound = boost::math::tools::root_epsilon<double>();
533 if (abs(x) >= taylor_2_bound) {
534 // approximation by taylor series in x at 0 up to order 4
535 result += x2 * x * (4.0 / 5.0);
536 }
537 }
538 return result;
539 }
540}

◆ calcPassiveMoveByJacobian() [1/2]

void calcPassiveMoveByJacobian ( const double &  byX,
const double &  byY,
TMatrixD &  jacobian,
const double  expandBelowChi = M_PI / 8 
) const

Calculate the jacobian matrix for the transport of the helix parameters, when moving the origin of the coordinate system to a new location.

Does not update the helix parameters in any way.

The jacobian matrix is written into the output parameter 'jacobian'. If output parameter is a 5x5 matrix only the derivatives of the 5 helix parameters are written If output parameter is a 6x6 matrix the derivatives of the 5 helix parameters and derivates of the two dimensional arc length are written. The derivatives of the arcLength2D are in the 6th row of the matrix.

Parameters
byXX displacement by which the origin of the coordinate system should be moved.
byYY displacement by which the origin of the coordinate system should be moved.
[out]jacobianThe jacobian matrix containing the derivatives of the five helix parameters after the move relative the orignal parameters.
expandBelowChiControl parameter below, which absolute value of chi an expansion of divergent terms shall be used. This parameter exists for testing the consistency of the expansion with the closed form. In other applications the parameter should remain at its default value.

Definition at line 309 of file Helix.cc.

313{
314 // 0. Preparations
315 // Initialise the return value to a unit matrix
316 jacobian.UnitMatrix();
317 assert(jacobian.GetNrows() == jacobian.GetNcols());
318 assert(jacobian.GetNrows() == 5 or jacobian.GetNrows() == 6);
319
320 // Fetch the helix parameters
321 const double omega = getOmega();
322 const double cosPhi0 = getCosPhi0();
323 const double sinPhi0 = getSinPhi0();
324 const double d0 = getD0();
325 const double tanLambda = getTanLambda();
326
327 // Prepare a delta vector, which is the vector from the perigee point to the new origin
328 // Split it in component parallel and a component orthogonal to tangent at the perigee.
329 const double deltaParallel = byX * cosPhi0 + byY * sinPhi0;
330 const double deltaOrthogonal = byY * cosPhi0 - byX * sinPhi0 + d0;
331 const double deltaCylindricalR = hypot(deltaOrthogonal, deltaParallel);
332 const double deltaCylindricalRSquared = deltaCylindricalR * deltaCylindricalR;
333
334 // Some commonly used terms - compare Karimaki 1990
335 const double A = 2 * deltaOrthogonal + omega * deltaCylindricalRSquared;
336 const double USquared = 1 + omega * A;
337 const double U = sqrt(USquared);
338 const double UOrthogonal = 1 + omega * deltaOrthogonal; // called nu in the Karimaki paper.
339 const double UParallel = omega * deltaParallel;
340 // Note U is a vector pointing from the middle of the projected circle scaled by a factor omega.
341 const double u = 1 + omega * d0; // just called u in the Karimaki paper.
342
343 // ---------------------------------------------------------------------------------------
344 // 1. Set the parts related to the xy coordinates
345 // Fills the upper left 3x3 corner of the jacobain
346
347 // a. Calculate the row related to omega
348 // The row related to omega is not different from the unit jacobian initialised above.
349 // jacobian(iOmega, iOmega) = 1.0; // From UnitMatrix above.
350
351 // b. Calculate the row related to d0
352 const double new_d0 = A / (1 + U);
353
354 jacobian(iD0, iOmega) = (deltaCylindricalR + new_d0) * (deltaCylindricalR - new_d0) / 2 / U;
355 jacobian(iD0, iPhi0) = - u * deltaParallel / U;
356 jacobian(iD0, iD0) = UOrthogonal / U;
357
358 // c. Calculate the row related to phi0
359 // Also calculate the derivatives of the arc length
360 // which are need for the row related to z.
361 const double dArcLength2D_dD0 = - UParallel / USquared;
362 const double dArcLength2D_dPhi0 = (omega * deltaCylindricalRSquared + deltaOrthogonal - d0 * UOrthogonal) / USquared;
363
364 jacobian(iPhi0, iD0) = - dArcLength2D_dD0 * omega;
365 jacobian(iPhi0, iPhi0) = u * UOrthogonal / USquared;
366 jacobian(iPhi0, iOmega) = -deltaParallel / USquared;
367
368 // For jacobian(iPhi0, iOmega) we have to dig deeper
369 // since the normal equations have a divergence for omega -> 0.
370 // This hinders not only the straight line case but extrapolations between points which are close together,
371 // like it happens when the current perigee is very close the new perigee.
372 // To have a smooth transition in this limit we have to carefully
373 // factor out the divergent quotients and approximate them with their taylor series.
374
375 const double chi = -atan2(UParallel, UOrthogonal);
376 double arcLength2D = 0;
377 double dArcLength2D_dOmega = 0;
378
379 if (fabs(chi) < std::min(expandBelowChi, M_PI / 2.0)) {
380 // Never expand for the far side of the circle by limiting the expandBelow to maximally half a pi.
381 // Close side of the circle
382 double principleArcLength2D = deltaParallel / UOrthogonal;
383 const double dPrincipleArcLength2D_dOmega = - principleArcLength2D * deltaOrthogonal / UOrthogonal;
384
385 const double x = principleArcLength2D * omega;
386 const double f = calcATanXDividedByX(x);
387 const double df_dx = calcDerivativeOfATanXDividedByX(x);
388
389 arcLength2D = principleArcLength2D * f;
390 dArcLength2D_dOmega = (f + x * df_dx) * dPrincipleArcLength2D_dOmega + principleArcLength2D * df_dx * principleArcLength2D;
391 } else {
392 // Far side of the circle
393 // If the far side of the circle is a well definied concept, omega is high enough that
394 // we can divide by it.
395 // Otherwise nothing can rescue us since the far side of the circle is so far away that no reasonable extrapolation can be made.
396 arcLength2D = - chi / omega;
397 dArcLength2D_dOmega = (-arcLength2D - jacobian(iPhi0, iOmega)) / omega;
398 }
399
400 // ---------------------------------------------------------------------------------------
401 // 2. Set the parts related to the z coordinate
402 // Since tanLambda stays the same there are no entries for tanLambda in the jacobian matrix
403 // For the new z0' = z0 + arcLength2D * tanLambda
404 jacobian(iZ0, iD0) = tanLambda * dArcLength2D_dD0;
405 jacobian(iZ0, iPhi0) = tanLambda * dArcLength2D_dPhi0;
406 jacobian(iZ0, iOmega) = tanLambda * dArcLength2D_dOmega;
407 jacobian(iZ0, iZ0) = 1.0; // From UnitMatrix above.
408 jacobian(iZ0, iTanLambda) = arcLength2D;
409
410 if (jacobian.GetNrows() == 6) {
411 // Also write the derivates of arcLength2D to the jacobian if the matrix size has space for them.
412 jacobian(iArcLength2D, iD0) = dArcLength2D_dD0;
413 jacobian(iArcLength2D, iPhi0) = dArcLength2D_dPhi0;
414 jacobian(iArcLength2D, iOmega) = dArcLength2D_dOmega;
415 jacobian(iArcLength2D, iZ0) = 0.0;
416 jacobian(iArcLength2D, iTanLambda) = 0;
417 jacobian(iArcLength2D, iArcLength2D) = 1.0;
418 }
419}
double getTanLambda() const
Getter for tan lambda, which is the z over two dimensional arc length slope of the track.
Definition: Helix.h:393
static double calcDerivativeOfATanXDividedByX(const double &x)
Implementation of the function d / dx (atan(x) / x) which handles small x values smoothly.
Definition: Helix.cc:509

◆ calcPassiveMoveByJacobian() [2/2]

TMatrixD calcPassiveMoveByJacobian ( const ROOT::Math::XYZVector &  by,
const double  expandBelowChi = M_PI / 8 
) const

Calculate the 5x5 jacobian matrix for the transport of the helix parameters, when moving the origin of the coordinate system to a new location.

Does not update the helix parameters in any way.

Parameters
byVector by which the origin of the coordinate system should be moved.
expandBelowChiControl parameter below, which absolute value of chi an expansion of divergent terms shall be used. This parameter exists for testing the consistency of the expansion with the closed form. In other applications the parameter should remain at its default value.
Returns
The jacobian matrix containing the derivatives of the five helix parameters after the move relative the orignal parameters.

Definition at line 301 of file Helix.cc.

302{
303 TMatrixD jacobian(5, 5);
304 calcPassiveMoveByJacobian(by.X(), by.Y(), jacobian, expandBelowChi);
305 return jacobian;
306}
TMatrixD calcPassiveMoveByJacobian(const ROOT::Math::XYZVector &by, const double expandBelowChi=M_PI/8) const
Calculate the 5x5 jacobian matrix for the transport of the helix parameters, when moving the origin o...
Definition: Helix.cc:301

◆ calcSecantLengthToArcLength2DFactor()

double calcSecantLengthToArcLength2DFactor ( const double &  secantLength2D) const

Helper function to calculate the factor between the dimensional secant length and the two dimensional arc length as seen in xy projection of the helix.

Definition at line 432 of file Helix.cc.

433{
434 double x = secantLength * getOmega() / 2.0;
435 return calcASinXDividedByX(x);
436}
static double calcASinXDividedByX(const double &x)
Implementation of the function asin(x) / x which handles small x values smoothly.
Definition: Helix.cc:438

◆ copyRelations()

void copyRelations ( const RelationsInterface< BASE > *  sourceObj)
inlineinherited

Copies all relations of sourceObj (pointing from or to sourceObj) to this object (including weights).

Useful if you want to make a complete copy of a StoreArray object to make modifications to it, but retain all information on linked objects.

Note: this only works if sourceObj inherits from the same base (e.g. RelationsObject), and only for related objects that also inherit from the same base.

Definition at line 170 of file RelationsObject.h.

171 {
172 if (!sourceObj)
173 return;
174 auto fromRels = sourceObj->getRelationsFrom<RelationsInterface<BASE>>("ALL");
175 for (unsigned int iRel = 0; iRel < fromRels.size(); iRel++) {
176 fromRels.object(iRel)->addRelationTo(this, fromRels.weight(iRel));
177 }
178
179 auto toRels = sourceObj->getRelationsTo<RelationsInterface<BASE>>("ALL");
180 for (unsigned int iRel = 0; iRel < toRels.size(); iRel++) {
181 this->addRelationTo(toRels.object(iRel), toRels.weight(iRel));
182 }
183 }
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).

◆ getAlpha()

double getAlpha ( const double  bZ)
static

Calculates the alpha value for a given magnetic field in Tesla.

Definition at line 109 of file Helix.cc.

110{
111 return 1.0 / (bZ * Const::speedOfLight) * 1E4;
112}
static const double speedOfLight
[cm/ns]
Definition: Const.h:695

◆ getArcLength2DAtCylindricalR()

double getArcLength2DAtCylindricalR ( const double &  cylindricalR) const

Calculates the transverse travel distance at the point the helix first reaches the given cylindrical radius.

Gives the two dimensional arc length in the forward direction that is traversed until a certain cylindrical radius is reached. Returns NAN, if the cylindrical radius can not be reached, either because it is to far outside or inside of the perigee.

Forward the result to getPositionAtArcLength2D() or getMomentumAtArcLength2D() in order to extrapolate to the cylinder detector boundaries.

The result always has a positive sign. Hence it refers to the forward direction. Adding a minus sign yields the point at the same cylindrical radius but in the backward direction.

Parameters
cylindricalRThe cylinder radius to extrapolate to.
Returns
The two dimensional arc length traversed to reach the cylindrical radius. NAN if it can not be reached.

Definition at line 128 of file Helix.cc.

129{
130 // Slight trick here
131 // Since the sought point is on the helix we treat it as the perigee
132 // and the origin as the point to extrapolate to.
133 // We know the distance of the origin to the circle, which is just d0
134 // The direct distance from the origin to the imaginary perigee is just the given cylindricalR.
135 const double dr = getD0();
136 const double deltaCylindricalR = cylindricalR;
137 const double absArcLength2D = calcArcLength2DAtDeltaCylindricalRAndDr(deltaCylindricalR, dr);
138 return absArcLength2D;
139}
double calcArcLength2DAtDeltaCylindricalRAndDr(const double &deltaCylindricalR, const double &dr) const
Helper method to calculate the two dimensional arc length from the perigee to a point at cylindrical ...
Definition: Helix.cc:581

◆ getArcLength2DAtNormalPlane()

double getArcLength2DAtNormalPlane ( const double &  x,
const double &  y,
const double &  nX,
const double &  nY 
) const

Calculates the arc length to reach a plane parallel to the z axes.

If there is no intersection with the plane NAN is returned

Parameters
xX coordinate of the support point of the plane
yY coordinate of the support point of the plane
nXX coordinate of the normal vector of the plane
nYY coordinate of the normal vector of the plane
Returns
Shortest two dimensional arc length to the plane

Definition at line 149 of file Helix.cc.

151{
152 // Construct the tangential vector to the plan in xy space
153 const double tX = nY;
154 const double tY = -nX;
155
156 // Fetch the helix parameters
157 const double omega = getOmega();
158 const double cosPhi0 = getCosPhi0();
159 const double sinPhi0 = getSinPhi0();
160 const double d0 = getD0();
161
162 // Prepare a delta vector, which is the vector from the perigee point to the support point of the plane
163 // Split it in component parallel and a component orthogonal to tangent at the perigee.
164 const double deltaParallel = byX * cosPhi0 + byY * sinPhi0;
165 const double deltaOrthogonal = byY * cosPhi0 - byX * sinPhi0 + d0;
166 const double deltaCylindricalR = hypot(deltaOrthogonal, deltaParallel);
167 const double deltaCylindricalRSquared = deltaCylindricalR * deltaCylindricalR;
168
169 const double UOrthogonal = 1 + omega * deltaOrthogonal; // called nu in the Karimaki paper.
170 const double UParallel = omega * deltaParallel;
171
172 // Some commonly used terms - compare Karimaki 1990
173 const double A = 2 * deltaOrthogonal + omega * deltaCylindricalRSquared;
174
175 const double tParallel = tX * cosPhi0 + tY * sinPhi0;
176 const double tOrthogonal = tY * cosPhi0 - tX * sinPhi0;
177 const double tCylindricalR = (tX * tX + tY * tY);
178
179 const double c = A / 2;
180 const double b = UParallel * tParallel + UOrthogonal * tOrthogonal;
181 const double a = omega / 2 * tCylindricalR;
182
183 const double discriminant = ((double)b) * b - 4 * a * c;
184 const double root = sqrt(discriminant);
185 const double bigSum = (b > 0) ? -b - root : -b + root;
186
187 const double bigOffset = bigSum / 2 / a;
188 const double smallOffset = 2 * c / bigSum;
189
190 const double distance1 = hypot(byX + bigOffset * tX, byY + bigOffset * tY);
191 const double distance2 = hypot(byX + smallOffset * tX, byY + smallOffset * tY);
192
193 if (distance1 < distance2) {
194 return getArcLength2DAtXY(byX + bigOffset * tX, byY + bigOffset * tY);
195 } else {
196 return getArcLength2DAtXY(byX + smallOffset * tX, byY + smallOffset * tY);
197 }
198}
const HepVector & a(void) const
Returns helix parameters.
Definition: Helix.h:428
double getArcLength2DAtXY(const double &x, const double &y) const
Calculates the two dimensional arc length at which the circle in the xy projection is closest to the ...
Definition: Helix.cc:141

◆ getArcLength2DAtXY()

double getArcLength2DAtXY ( const double &  x,
const double &  y 
) const

Calculates the two dimensional arc length at which the circle in the xy projection is closest to the point.

This calculates the dimensional arc length to the closest approach in xy projection. Hence, it only optimizes the distance in x and y. This is sufficent to extrapolate to an axial wire. For stereo wires this is not optimal, since the distance in the z direction also plays a role.

Parameters
xX coordinate of the point to which to extrapolate
yY coordinate of the point to which to extrapolate
Returns
The two dimensional arc length from the perigee at which the closest approach is reached

Definition at line 141 of file Helix.cc.

142{
143 double dr = 0;
144 double arcLength2D = 0;
145 calcArcLength2DAndDrAtXY(x, y, arcLength2D, dr);
146 return arcLength2D;
147}
void calcArcLength2DAndDrAtXY(const double &x, const double &y, double &arcLength2D, double &dr) const
Helper method to calculate the signed two dimensional arc length and the signed distance to the circl...
Definition: Helix.cc:545

◆ getArrayIndex()

int getArrayIndex ( ) const
inlineinherited

Returns this object's array index (in StoreArray), or -1 if not found.

Definition at line 385 of file RelationsObject.h.

386 {
388 return m_cacheArrayIndex;
389 }
bool findStoreEntry(const TObject *object, StoreEntry *&entry, int &index)
Find an object in an array in the data store.
Definition: DataStore.cc:398

◆ getArrayName()

std::string getArrayName ( ) const
inlineinherited

Get name of array this object is stored in, or "" if not found.

Definition at line 377 of file RelationsObject.h.

◆ getArrayPointer()

TClonesArray * getArrayPointer ( ) const
inlineprotectedinherited

Returns the pointer to the raw DataStore array holding this object (protected since these arrays are easy to misuse).

Definition at line 418 of file RelationsObject.h.

419 {
422 return nullptr;
424 }
TClonesArray * getPtrAsArray() const
Return ptr cast to TClonesArray.
Definition: StoreEntry.cc:83

◆ getChargeSign()

short getChargeSign ( ) const

Return track charge sign (1, 0 or -1).

Definition at line 114 of file Helix.cc.

115{
116 return boost::math::sign(getOmega());
117}

◆ getCosPhi0()

double getCosPhi0 ( ) const
inline

Getter for the cosine of the azimuth angle of travel direction at the perigee.

Definition at line 381 of file Helix.h.

381{ return std::cos(getPhi0()); }
double getPhi0() const
Getter for phi0, which is the azimuth angle of the transverse momentum at the perigee.
Definition: Helix.h:378

◆ getCotTheta()

double getCotTheta ( ) const
inline

Getter for cot theta, which is the z over two dimensional arc length slope of the track.

Synomym to tan lambda.

Definition at line 396 of file Helix.h.

396{ return m_tanLambda; }

◆ getD0()

double getD0 ( ) const
inline

Getter for d0, which is the signed distance to the perigee in the r-phi plane.

The signed distance from the origin to the perigee. The sign is positive (negative), if the angle from the xy perigee position vector to the transverse momentum vector is +pi/2 (-pi/2). d0 has the same sign as getPerigee().Cross(getMomentum()).Z().

Definition at line 372 of file Helix.h.

372{ return m_d0; }

◆ getDirection()

ROOT::Math::XYZVector getDirection ( ) const

Getter for unit vector of momentum at the perigee position.

This is mainly useful cases where curvature is zero (pT is infinite)

Definition at line 94 of file Helix.cc.

95{
96 return ROOT::Math::XYZVector(getCosPhi0(), getSinPhi0(), getTanLambda());
97}

◆ getInfo()

std::string getInfo ( ) const
inlineinherited

Return a short summary of this object's contents in raw text format.

Returns the contents of getInfoHTML() while translating line-breaks etc.

Note
: You don't need to implement this function (it's not virtual), getInfoHTML() is enough.

Definition at line 370 of file RelationsObject.h.

371 {
373 }
virtual std::string getInfoHTML() const
Return a short summary of this object's contents in HTML format.
std::string htmlToPlainText(const std::string &html)
See RelationsObject::getInfo()

◆ getInfoHTML()

virtual std::string getInfoHTML ( ) const
inlinevirtualinherited

Return a short summary of this object's contents in HTML format.

Reimplement this in your own class to provide useful output for display or debugging purposes. For example, you might do something like:

std::stringstream out;
out << "<b>PDG</b>: " << m_pdg << "<br>";
out << "<b>Covariance Matrix</b>: " << HTML::getString(getCovariance5()) << "<br>";
return out.str();
std::string getString(const TMatrixFBase &matrix, int precision=2, bool color=true)
get HTML table representing a matrix.
Definition: HTML.cc:24
See also
Particle::getInfoHTML() for a more complex example.
HTML for some utility functions.
Use getInfo() to get a raw text version of this output.

Reimplemented in Particle, Cluster, MCParticle, PIDLikelihood, SoftwareTriggerResult, Track, TrackFitResult, TRGSummary, and RecoTrack.

Definition at line 362 of file RelationsObject.h.

362{ return ""; }

◆ getKappa()

double getKappa ( const double  bZ) const

Getter for kappa, which is charge / transverse momentum or equivalently omega * alpha.

Definition at line 104 of file Helix.cc.

105{
106 return getOmega() * getAlpha(bZ);
107}
static double getAlpha(const double bZ)
Calculates the alpha value for a given magnetic field in Tesla.
Definition: Helix.cc:109

◆ getMomentum()

ROOT::Math::XYZVector getMomentum ( const double  bZ) const

Getter for vector of momentum at the perigee position.

As we calculate recalculate the momentum from a geometric helix, we need an estimate of the magnetic field along the z-axis to give back the momentum.

Parameters
bZMagnetic field at the perigee.

Definition at line 89 of file Helix.cc.

90{
91 return ROOT::Math::XYZVector(getMomentumX(bZ), getMomentumY(bZ), getMomentumZ(bZ));
92}
double getMomentumZ(const double bZ) const
Calculates the z momentum of the particle at the perigee point.
Definition: Helix.cc:84
double getMomentumY(const double bZ) const
Calculates the y momentum of the particle at the perigee point.
Definition: Helix.cc:79
double getMomentumX(const double bZ) const
Calculates the x momentum of the particle at the perigee point.
Definition: Helix.cc:74

◆ getMomentumAtArcLength2D()

ROOT::Math::XYZVector getMomentumAtArcLength2D ( const double &  arcLength2D,
const double &  bz 
) const

Calculates the momentum vector at the given two dimensional arc length.

Parameters
arcLength2DTwo dimensional arc length to be traversed.
bzMagnetic field strength in the z direction.

Definition at line 270 of file Helix.cc.

271{
272 ROOT::Math::XYZVector momentum = getTangentialAtArcLength2D(arcLength2D);
273 const double pr = getTransverseMomentum(bz);
274 momentum *= pr;
275
276 return momentum;
277}
ROOT::Math::XYZVector getTangentialAtArcLength2D(const double &arcLength2D) const
Calculates the tangential vector to the helix curve at the given two dimensional arc length.
Definition: Helix.cc:245
double getTransverseMomentum(const double bZ) const
Getter for the absolute value of the transverse momentum at the perigee.
Definition: Helix.cc:99

◆ getMomentumX()

double getMomentumX ( const double  bZ) const

Calculates the x momentum of the particle at the perigee point.

Parameters
bZZ component of the magnetic field in Tesla

Definition at line 74 of file Helix.cc.

75{
76 return getCosPhi0() * getTransverseMomentum(bZ);
77}

◆ getMomentumY()

double getMomentumY ( const double  bZ) const

Calculates the y momentum of the particle at the perigee point.

Parameters
bZZ component of the magnetic field in Tesla

Definition at line 79 of file Helix.cc.

80{
81 return getSinPhi0() * getTransverseMomentum(bZ);
82}

◆ getMomentumZ()

double getMomentumZ ( const double  bZ) const

Calculates the z momentum of the particle at the perigee point.

Parameters
bZZ component of the magnetic field in Tesla

Definition at line 84 of file Helix.cc.

85{
87}

◆ getName()

virtual std::string getName ( ) const
inlinevirtualinherited

Return a short name that describes this object, e.g.

pi+ for an MCParticle.

Reimplemented in Particle, MCParticle, and SpacePoint.

Definition at line 344 of file RelationsObject.h.

344{ return ""; }

◆ getOmega()

double getOmega ( ) const
inline

Getter for omega, which is a signed curvature measure of the track.

The sign is equivalent to the charge of the particle.

Definition at line 387 of file Helix.h.

387{ return m_omega; }

◆ getPerigee()

ROOT::Math::XYZVector getPerigee ( ) const

Getter for the perigee position.

Definition at line 69 of file Helix.cc.

70{
71 return ROOT::Math::XYZVector(getPerigeeX(), getPerigeeY(), getPerigeeZ());
72}
double getPerigeeX() const
Calculates the x coordinate of the perigee point.
Definition: Helix.cc:54
double getPerigeeZ() const
Calculates the z coordinate of the perigee point.
Definition: Helix.cc:64
double getPerigeeY() const
Calculates the y coordinate of the perigee point.
Definition: Helix.cc:59

◆ getPerigeeX()

double getPerigeeX ( ) const

Calculates the x coordinate of the perigee point.

Definition at line 54 of file Helix.cc.

55{
56 return getD0() * getSinPhi0();
57}

◆ getPerigeeY()

double getPerigeeY ( ) const

Calculates the y coordinate of the perigee point.

Definition at line 59 of file Helix.cc.

60{
61 return -getD0() * getCosPhi0();
62}

◆ getPerigeeZ()

double getPerigeeZ ( ) const

Calculates the z coordinate of the perigee point.

Definition at line 64 of file Helix.cc.

65{
66 return getZ0();
67}
double getZ0() const
Getter for z0, which is the z coordinate of the perigee.
Definition: Helix.h:390

◆ getPhi0()

double getPhi0 ( ) const
inline

Getter for phi0, which is the azimuth angle of the transverse momentum at the perigee.

getMomentum().Phi() == getPhi0() holds.

Definition at line 378 of file Helix.h.

378{ return m_phi0; }

◆ getPositionAtArcLength2D()

ROOT::Math::XYZVector getPositionAtArcLength2D ( const double &  arcLength2D) const

Calculates the position on the helix at the given two dimensional arc length.

Parameters
arcLength2DTwo dimensional arc length to be traversed.

Definition at line 201 of file Helix.cc.

202{
203 /*
204 / \ / \ / \
205 | x | | cos phi0 -sin phi0 | | - sin(chi) / omega |
206 | | = | | * | |
207 | y | | sin phi0 cos phi0 | | -(1 - cos(chi)) / omega - d0 |
208 \ / \ / \ /
209
210 and
211
212 z = tanLambda * arclength + z0;
213
214 where chi = -arcLength2D * omega
215
216 Here arcLength2D means the arc length of the circle in the xy projection
217 traversed in the forward direction.
218 */
219
220 // First calculating the position assuming the circle center lies on the y axes (phi0 = 0)
221 // Rotate to the right phi position afterwards
222 // Using the sinus cardinalis yields expressions that are smooth in the limit of omega -> 0
223
224 // Do calculations in double
225 const double chi = -arcLength2D * getOmega();
226 const double chiHalf = chi / 2.0;
227
228 using boost::math::sinc_pi;
229
230 const double x = arcLength2D * sinc_pi(chi);
231 const double y = arcLength2D * sinc_pi(chiHalf) * sin(chiHalf) - getD0();
232
233 // const double z = s * tan lambda + z0
234 const double z = fma((double)arcLength2D, getTanLambda(), getZ0());
235
236 // Unrotated position
237 ROOT::Math::XYZVector position(x, y, z);
238
239 // Rotate to the right phi0 position
240 position = ROOT::Math::VectorUtil::RotateZ(position, getPhi0());
241
242 return position;
243}

◆ getRelated()

T * getRelated ( const std::string &  name = "",
const std::string &  namedRelation = "" 
) const
inlineinherited

Get the object to or from which this object has a relation.

Template Parameters
TThe class of objects to or from which the relation points.
Parameters
nameThe name of the store array to or from which the relation points. If empty the default store array name for class T will be used. If the special name "ALL" is given all store arrays containing objects of type T are considered.
namedRelationAdditional name for the relation, or "" for the default naming
Returns
The first related object or a null pointer.

Definition at line 278 of file RelationsObject.h.

279 {
281 T::Class(), name, namedRelation).object);
282 }
@ c_BothSides
Combination of c_FromSide and c_ToSide.
Definition: DataStore.h:79
Belle2::RelationEntry getRelationWith(ESearchSide searchSide, const TObject *object, StoreEntry *&entry, int &index, const TClass *withClass, const std::string &withName, const std::string &namedRelation)
Get the first relation between an object and another object in a store array.
Definition: DataStore.cc:597
TObject * object
Pointer to the object.
Definition: RelationEntry.h:32

◆ getRelatedFrom()

FROM * getRelatedFrom ( const std::string &  name = "",
const std::string &  namedRelation = "" 
) const
inlineinherited

Get the object from which this object has a relation.

Template Parameters
FROMThe class of objects from which the relation points.
Parameters
nameThe name of the store array from which the relation points. If empty the default store array name for class FROM will be used. If the special name "ALL" is given all store arrays containing objects of type FROM are considered.
namedRelationAdditional name for the relation, or "" for the default naming
Returns
The first related object or a null pointer.

Definition at line 263 of file RelationsObject.h.

264 {
266 m_cacheArrayIndex, FROM::Class(), name, namedRelation).object);
267 }
@ c_FromSide
Return relations/objects pointed from (to a given object).
Definition: DataStore.h:77

◆ getRelatedFromWithWeight()

std::pair< FROM *, float > getRelatedFromWithWeight ( const std::string &  name = "",
const std::string &  namedRelation = "" 
) const
inlineinherited

Get first related object & weight of relation pointing from an array.

Template Parameters
FROMThe class of objects from which the relation points.
Parameters
nameThe name of the store array from which the relation points. If empty the default store array name for class FROM will be used. If the special name "ALL" is given all store arrays containing objects of type FROM are considered.
namedRelationAdditional name for the relation, or "" for the default naming
Returns
Pair of first related object and the relation weight, or (NULL, 1.0) if none found.

Definition at line 314 of file RelationsObject.h.

316 {
318 FROM::Class(), name, namedRelation);
319 return std::make_pair(static_cast<FROM*>(entry.object), entry.weight);
320 }

◆ getRelatedTo()

TO * getRelatedTo ( const std::string &  name = "",
const std::string &  namedRelation = "" 
) const
inlineinherited

Get the object to which this object has a relation.

Template Parameters
TOThe class of objects to which the relation points.
Parameters
nameThe name of the store array to which the relation points. If empty the default store array name for class TO will be used. If the special name "ALL" is given all store arrays containing objects of type TO are considered.
namedRelationAdditional name for the relation, or "" for the default naming
Returns
The first related object or a null pointer.

Definition at line 248 of file RelationsObject.h.

249 {
251 TO::Class(), name, namedRelation).object);
252 }
@ c_ToSide
Return relations/objects pointed to (from a given object).
Definition: DataStore.h:78

◆ getRelatedToWithWeight()

std::pair< TO *, float > getRelatedToWithWeight ( const std::string &  name = "",
const std::string &  namedRelation = "" 
) const
inlineinherited

Get first related object & weight of relation pointing to an array.

Template Parameters
TOThe class of objects to which the relation points.
Parameters
nameThe name of the store array to which the relation points. If empty the default store array name for class TO will be used. If the special name "ALL" is given all store arrays containing objects of type TO are considered.
namedRelationAdditional name for the relation, or "" for the default naming
Returns
Pair of first related object and the relation weight, or (NULL, 1.0) if none found.

Definition at line 297 of file RelationsObject.h.

299 {
301 TO::Class(), name, namedRelation);
302 return std::make_pair(static_cast<TO*>(entry.object), entry.weight);
303 }

◆ getRelatedWithWeight()

std::pair< T *, float > getRelatedWithWeight ( const std::string &  name = "",
const std::string &  namedRelation = "" 
) const
inlineinherited

Get first related object & weight of relation pointing from/to an array.

Template Parameters
TThe class of objects to or from which the relation points.
Parameters
nameThe name of the store array to or from which the relation points. If empty the default store array name for class T will be used. If the special name "ALL" is given all store arrays containing objects of type T are considered.
namedRelationAdditional name for the relation, or "" for the default naming
Returns
Pair of first related object and the relation weight, or (NULL, 1.0) if none found.

Definition at line 331 of file RelationsObject.h.

333 {
335 T::Class(), name, namedRelation);
336 return std::make_pair(static_cast<T*>(entry.object), entry.weight);
337 }

◆ getRelationsFrom()

RelationVector< FROM > getRelationsFrom ( const std::string &  name = "",
const std::string &  namedRelation = "" 
) const
inlineinherited

Get the relations that point from another store array to this object.

Template Parameters
FROMThe class of objects from which the relations point.
Parameters
nameThe name of the store array from which the relations point. If empty the default store array name for class FROM will be used. If the special name "ALL" is given all store arrays containing objects of type FROM are considered.
namedRelationAdditional name for the relation, or "" for the default naming
Returns
A vector of relations.

Definition at line 212 of file RelationsObject.h.

214 {
216 m_cacheArrayIndex, FROM::Class(), name, namedRelation));
217 }
RelationVector< T > getRelationsWith(const std::string &name="", const std::string &namedRelation="") const
Get the relations between this object and another store array.

◆ getRelationsTo()

RelationVector< TO > getRelationsTo ( const std::string &  name = "",
const std::string &  namedRelation = "" 
) const
inlineinherited

Get the relations that point from this object to another store array.

Template Parameters
TOThe class of objects to which the relations point.
Parameters
nameThe name of the store array to which the relations point. If empty the default store array name for class TO will be used. If the special name "ALL" is given all store arrays containing objects of type TO are considered.
namedRelationAdditional name for the relation, or "" for the default naming
Returns
A vector of relations.

Definition at line 197 of file RelationsObject.h.

198 {
200 m_cacheArrayIndex, TO::Class(), name, namedRelation));
201 }

◆ getRelationsWith()

RelationVector< T > getRelationsWith ( const std::string &  name = "",
const std::string &  namedRelation = "" 
) const
inlineinherited

Get the relations between this object and another store array.

Relations in both directions are returned.

Template Parameters
TThe class of objects to or from which the relations point.
Parameters
nameThe name of the store array to or from which the relations point. If empty the default store array name for class T will be used. If the special name "ALL" is given all store arrays containing objects of type T are considered.
namedRelationAdditional name for the relation, or "" for the default naming
Returns
A vector of relations.

Definition at line 230 of file RelationsObject.h.

231 {
233 m_cacheArrayIndex, T::Class(), name, namedRelation));
234 }

◆ getSinPhi0()

double getSinPhi0 ( ) const
inline

Getter for the cosine of the azimuth angle of travel direction at the perigee.

Definition at line 384 of file Helix.h.

384{ return std::sin(getPhi0()); }

◆ getTangentialAtArcLength2D()

ROOT::Math::XYZVector getTangentialAtArcLength2D ( const double &  arcLength2D) const

Calculates the tangential vector to the helix curve at the given two dimensional arc length.

The tangential vector is the derivative of the position with respect to the two dimensional arc length It is normalised such that the cylindrical radius of the result is 1

getTangentialAtArcLength2D(arcLength2D).Perp() == 1 holds.

Parameters
arcLength2DTwo dimensional arc length to be traversed.
Returns
Tangential vector normalised to unit transverse component / cylindrical radius.

Definition at line 245 of file Helix.cc.

246{
247 const double omega = getOmega();
248 const double phi0 = getPhi0();
249 const double tanLambda = getTanLambda();
250
251 const double chi = - omega * arcLength2D;
252
253 const double tx = cos(chi + phi0);
254 const double ty = sin(chi + phi0);
255 const double tz = tanLambda;
256
257 ROOT::Math::XYZVector tangential(tx, ty, tz);
258 return tangential;
259}

◆ getTanLambda()

double getTanLambda ( ) const
inline

Getter for tan lambda, which is the z over two dimensional arc length slope of the track.

Definition at line 393 of file Helix.h.

393{ return m_tanLambda; }

◆ getTransverseMomentum()

double getTransverseMomentum ( const double  bZ) const

Getter for the absolute value of the transverse momentum at the perigee.

Parameters
bZMagnetic field at the perigee

Definition at line 99 of file Helix.cc.

100{
101 return 1 / std::fabs(getAlpha(bZ) * getOmega());
102}

◆ getUnitTangentialAtArcLength2D()

ROOT::Math::XYZVector getUnitTangentialAtArcLength2D ( const double &  arcLength2D) const

Calculates the unit tangential vector to the helix curve at the given two dimensional arc length.

Parameters
arcLength2DTwo dimensional arc length to be traversed.

Definition at line 261 of file Helix.cc.

262{
263 ROOT::Math::XYZVector unitTangential = getTangentialAtArcLength2D(arcLength2D);
264 const double norm = hypot(1, getTanLambda());
265 const double invNorm = 1 / norm;
266 unitTangential *= invNorm;
267 return unitTangential;
268}

◆ getZ0()

double getZ0 ( ) const
inline

Getter for z0, which is the z coordinate of the perigee.

Definition at line 390 of file Helix.h.

390{ return m_z0; }

◆ passiveMoveBy() [1/2]

double passiveMoveBy ( const double &  byX,
const double &  byY,
const double &  byZ 
)

Moves origin of the coordinate system (passive transformation) orthogonal to the z axis by the given vector.

Updates the helix inplace.

Parameters
byXX displacement by which the origin of the coordinate system should be moved.
byYY displacement by which the origin of the coordinate system should be moved.
byZZ displacement by which the origin of the coordinate system should be moved.
Returns
The double value is the two dimensional arc length, which has the be traversed from the old perigee to the new.

Update the parameters inplace. Omega and tan lambda are unchanged

Definition at line 279 of file Helix.cc.

282{
283 // First calculate the distance of the new origin to the helix in the xy projection
284 double new_d0 = 0;
285 double arcLength2D = 0;
286 calcArcLength2DAndDrAtXY(byX, byY, arcLength2D, new_d0);
287
288 // Third the new phi0 and z0 can be calculated from the arc length
289 double chi = -arcLength2D * getOmega();
290 double new_phi0 = m_phi0 + chi;
291 double new_z0 = getZ0() - byZ + getTanLambda() * arcLength2D;
292
294 m_d0 = new_d0;
295 m_phi0 = new_phi0;
296 m_z0 = new_z0;
297
298 return arcLength2D;
299}

◆ passiveMoveBy() [2/2]

double passiveMoveBy ( const ROOT::Math::XYZVector &  by)
inline

Moves origin of the coordinate system (passive transformation) by the given vector.

Updates the helix inplace.

Parameters
byVector by which the origin of the coordinate system should be moved.
Returns
The double value is the two dimensional arc length, which has the be traversed from the old perigee to the new.

Definition at line 245 of file Helix.h.

246 { return passiveMoveBy(by.X(), by.Y(), by.Z()); }
double passiveMoveBy(const ROOT::Math::XYZVector &by)
Moves origin of the coordinate system (passive transformation) by the given vector.
Definition: Helix.h:245

◆ reverse()

void reverse ( )

Reverses the direction of travel of the helix in place.

The same points that are located on the helix stay are the same after the transformation, but have the opposite two dimensional arc length associated to them. The momentum at each point is reversed. The charge sign is changed to its opposite by this transformation.

◆ reversePhi()

double reversePhi ( const double &  phi)
static

Reverses an azimuthal angle to the opposite direction.

Parameters
phiA angle in [-pi, pi]
Returns
The angle for the opposite direction in [-pi, pi]

Definition at line 421 of file Helix.cc.

422{
423 return std::remainder(phi + M_PI, 2 * M_PI);
424}

◆ setCartesian()

void setCartesian ( const ROOT::Math::XYZVector &  position,
const ROOT::Math::XYZVector &  momentum,
const short int  charge,
const double  bZ 
)
private

Cartesian to Perigee conversion.


Definition at line 588 of file Helix.cc.

592{
593 assert(abs(charge) <= 1); // Not prepared for doubly-charged particles.
594 const double alpha = getAlpha(bZ);
595
596 // We allow for the case that position, momentum are not given
597 // exactly in the perigee. Therefore we have to transform the momentum
598 // with the position as the reference point and then move the coordinate system
599 // to the origin.
600
601 const double x = position.X();
602 const double y = position.Y();
603 const double z = position.Z();
604
605 const double px = momentum.X();
606 const double py = momentum.Y();
607 const double pz = momentum.Z();
608
609 const double ptinv = 1 / hypot(px, py);
610 const double omega = charge * ptinv / alpha;
611 const double tanLambda = ptinv * pz;
612 const double phi0 = atan2(py, px);
613 const double z0 = z;
614 const double d0 = 0;
615
616 m_omega = omega;
617 m_phi0 = phi0;
618 m_d0 = d0;
619 m_tanLambda = tanLambda;
620 m_z0 = z0;
621
622 passiveMoveBy(-x, -y, 0);
623}
double charge(int pdgCode)
Returns electric charge of a particle with given pdg code.
Definition: EvtPDLUtil.cc:44

Friends And Related Function Documentation

◆ operator<<

std::ostream & operator<< ( std::ostream &  output,
const Helix helix 
)
friend

Output operator for debugging and the generation of unittest error messages.

Definition at line 631 of file Helix.cc.

632 {
633 return output
634 << "Helix("
635 << "d0=" << helix.getD0() << ", "
636 << "phi0=" << helix.getPhi0() << ", "
637 << "omega=" << helix.getOmega() << ", "
638 << "z0=" << helix.getZ0() << ", "
639 << "tanLambda=" << helix.getTanLambda() << ")";
640 }

Member Data Documentation

◆ m_cacheArrayIndex

int m_cacheArrayIndex
mutableprivateinherited

Cache of the index in the TClonesArray to which this object belongs.

Definition at line 432 of file RelationsObject.h.

◆ m_cacheDataStoreEntry

DataStore::StoreEntry* m_cacheDataStoreEntry
mutableprivateinherited

Cache of the data store entry to which this object belongs.

Definition at line 429 of file RelationsObject.h.

◆ m_d0

Double32_t m_d0
private

Memory for the signed distance to the perigee.

The sign is the same as of the z component of getPerigee().Cross(getMomentum()).

Definition at line 414 of file Helix.h.

◆ m_omega

Double32_t m_omega
private

Memory for the curvature of the signed curvature.

Definition at line 420 of file Helix.h.

◆ m_phi0

Double32_t m_phi0
private

Memory for the azimuth angle between the transverse momentum and the x axis, which is in [-pi, pi].

Definition at line 417 of file Helix.h.

◆ m_tanLambda

Double32_t m_tanLambda
private

Memory for the slope of the track in the z coordinate over the two dimensional arc length (dz/ds)

Definition at line 426 of file Helix.h.

◆ m_z0

Double32_t m_z0
private

Memory for the z coordinate of the perigee.

Definition at line 423 of file Helix.h.


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