ACTS
Experiment-independent tracking
Loading...
Searching...
No Matches
Acts::VectorHelpers Namespace Reference

Functions

double cast (const Vector3 &position, AxisDirection aDir)
 Helper method to extract the binning value from a 3D vector.
SquareMatrix3 cross (const SquareMatrix3 &m, const Vector3 &v)
 Calculates column-wise cross products of a matrix and a vector and stores the result column-wise in a matrix.
template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime == 3)
double deltaR (const Eigen::MatrixBase< Derived > &v1, const Eigen::MatrixBase< Derived > &v2)
 Calculate the deltaR between two vectors.
template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime == -1)
double eta (const Eigen::MatrixBase< Derived > &v) noexcept
 Calculate the pseudorapidity for a vector with dynamic size.
template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime == 3)
double eta (const Eigen::MatrixBase< Derived > &v) noexcept
 Calculate the pseudorapidity for a vector with static size.
std::array< double, 4 > evaluateTrigonomics (const Vector3 &direction)
 Fast evaluation of trigonomic functions.
std::pair< double, double > incidentAngles (const Acts::Vector3 &direction, const Acts::RotationMatrix3 &globalToLocal)
 Calculate the incident angles of a vector with in a given reference frame.
template<typename vector3_t>
auto makeVector4 (const Eigen::MatrixBase< vector3_t > &vec3, typename vector3_t::Scalar w) -> Eigen::Matrix< typename vector3_t::Scalar, 4, 1 >
 Construct a four-vector from a three-vector and scalar fourth component.
template<typename Derived>
double perp (const Eigen::MatrixBase< Derived > &v) noexcept
 Calculate radius in the transverse (xy) plane of a vector.
template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime == -1)
double phi (const Eigen::MatrixBase< Derived > &v) noexcept
 Calculate phi (transverse plane angle) from compatible Eigen types with dynamic size.
template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime >= 2)
double phi (const Eigen::MatrixBase< Derived > &v) noexcept
 Calculate phi (transverse plane angle) from compatible Eigen types with static size.
template<typename T>
requires requires { { v.phi() } -> std::floating_point; }
double phi (const T &v) noexcept
 Calculate phi (transverse plane angle) from anything implementing a method like phi() returning anything convertible to double.
auto position (const FreeVector &params)
 Access the three-position components in a free parameters vector.
auto position (const Vector4 &pos4)
 Access the three-position components in a four-position vector.
template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime == -1)
double theta (const Eigen::MatrixBase< Derived > &v) noexcept
 Calculate the theta angle (longitudinal w.r.t.
template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime == 3)
double theta (const Eigen::MatrixBase< Derived > &v) noexcept
 Calculate the theta angle (longitudinal w.r.t.

Function Documentation

◆ cast()

double Acts::VectorHelpers::cast ( const Vector3 & position,
AxisDirection aDir )

Helper method to extract the binning value from a 3D vector.

For this method a 3D vector is required to guarantee all potential axis directions to be casted from

Parameters
positionis the position in global
aDiris the axis direction to be extracted
Returns
the value of the binning direction

◆ cross()

SquareMatrix3 Acts::VectorHelpers::cross ( const SquareMatrix3 & m,
const Vector3 & v )

Calculates column-wise cross products of a matrix and a vector and stores the result column-wise in a matrix.

Parameters
[in]mMatrix that will be used for cross products
[in]vVector for cross products
Returns
Constructed matrix

◆ deltaR()

template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime == 3)
double Acts::VectorHelpers::deltaR ( const Eigen::MatrixBase< Derived > & v1,
const Eigen::MatrixBase< Derived > & v2 )

Calculate the deltaR between two vectors.

Note
DeltaR is defined as sqrt(deltaPhi^2 + deltaEta^2)
Template Parameters
DerivedEigen derived concrete type
Parameters
v1First vector
v2Second vector
Returns
The deltaR value

◆ eta() [1/2]

template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime == -1)
double Acts::VectorHelpers::eta ( const Eigen::MatrixBase< Derived > & v)
noexcept

Calculate the pseudorapidity for a vector with dynamic size.

Template Parameters
DerivedEigen derived concrete type with dynamic size
Parameters
vAny vector like Eigen type with dynamic size
Note
Will abort execution if the number of rows of v is not exactly 3.
Returns
The pseudorapidity value

◆ eta() [2/2]

template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime == 3)
double Acts::VectorHelpers::eta ( const Eigen::MatrixBase< Derived > & v)
noexcept

Calculate the pseudorapidity for a vector with static size.

Template Parameters
DerivedEigen derived concrete type with compile-time size == 3
Parameters
vAny 3D vector like Eigen type with static size
Returns
The pseudorapidity value

◆ evaluateTrigonomics()

std::array< double, 4 > Acts::VectorHelpers::evaluateTrigonomics ( const Vector3 & direction)

Fast evaluation of trigonomic functions.

Parameters
directionfor this evaluatoin
Returns
cos(phi), sin(phi), cos(theta), sin(theta), 1/sin(theta)

◆ incidentAngles()

std::pair< double, double > Acts::VectorHelpers::incidentAngles ( const Acts::Vector3 & direction,
const Acts::RotationMatrix3 & globalToLocal )

Calculate the incident angles of a vector with in a given reference frame.

Template Parameters
DerivedEigen derived concrete type
Parameters
directionThe crossing direction in the global frame
globalToLocalRotation from global to local frame
Returns
The angles of incidence in the two normal planes

◆ makeVector4()

template<typename vector3_t>
auto Acts::VectorHelpers::makeVector4 ( const Eigen::MatrixBase< vector3_t > & vec3,
typename vector3_t::Scalar w )->Eigen::Matrix< typenamevector3_t::Scalar, 4, 1 >

Construct a four-vector from a three-vector and scalar fourth component.

◆ perp()

template<typename Derived>
double Acts::VectorHelpers::perp ( const Eigen::MatrixBase< Derived > & v)
noexcept

Calculate radius in the transverse (xy) plane of a vector.

Template Parameters
DerivedEigen derived concrete type
Parameters
vAny vector like Eigen type, static or dynamic
Note
Will static assert that the number of rows of v is at least 2, or in case of dynamic size, will abort execution if that is not the case.
Returns
The transverse radius value.

◆ phi() [1/3]

template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime == -1)
double Acts::VectorHelpers::phi ( const Eigen::MatrixBase< Derived > & v)
noexcept

Calculate phi (transverse plane angle) from compatible Eigen types with dynamic size.

Template Parameters
DerivedEigen derived concrete type with dynamic size
Parameters
vAny vector like Eigen type with dynamic size
Note
Will abort execution if the number of rows of v is less than 2.
Returns
The value of the angle in the transverse plane.

◆ phi() [2/3]

template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime >= 2)
double Acts::VectorHelpers::phi ( const Eigen::MatrixBase< Derived > & v)
noexcept

Calculate phi (transverse plane angle) from compatible Eigen types with static size.

Template Parameters
DerivedEigen derived concrete type with compile-time known size >= 2
Parameters
vAny vector like Eigen type with static size
Returns
The value of the angle in the transverse plane.

◆ phi() [3/3]

template<typename T>
requires requires { { v.phi() } -> std::floating_point; }
double Acts::VectorHelpers::phi ( const T & v)
noexcept

Calculate phi (transverse plane angle) from anything implementing a method like phi() returning anything convertible to double.

Template Parameters
Tanything that has a phi method
Parameters
vAny type that implements a phi method
Returns
The phi value

◆ position() [1/2]

auto Acts::VectorHelpers::position ( const FreeVector & params)

Access the three-position components in a free parameters vector.

◆ position() [2/2]

auto Acts::VectorHelpers::position ( const Vector4 & pos4)

Access the three-position components in a four-position vector.

◆ theta() [1/2]

template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime == -1)
double Acts::VectorHelpers::theta ( const Eigen::MatrixBase< Derived > & v)
noexcept

Calculate the theta angle (longitudinal w.r.t.

z axis) of a vector with dynamic size

Template Parameters
DerivedEigen derived concrete type with dynamic size
Parameters
vAny vector like Eigen type with dynamic size
Note
Will abort execution if the number of rows of v is not exactly 3.
Returns
The theta value

◆ theta() [2/2]

template<typename Derived>
requires (Eigen::MatrixBase<Derived>::RowsAtCompileTime == 3)
double Acts::VectorHelpers::theta ( const Eigen::MatrixBase< Derived > & v)
noexcept

Calculate the theta angle (longitudinal w.r.t.

z axis) of a vector with static size

Template Parameters
DerivedEigen derived concrete type with compile-time size == 3
Parameters
vAny 3D vector like Eigen type with static size
Returns
The theta value