math Namespace Reference
Kynema API
A flexible multibody structural dynamics code for wind turbines
|
Functions | |
std::vector< double > | GetGllLocations (size_t order) |
std::vector< double > | GetGllWeights (size_t order) |
void | LinearInterpWeights (double x, std::span< const double > xs, std::vector< double > &weights) |
Computes weights for linear interpolation. | |
double | LinearInterp (double x, std::span< const double > xs, std::span< const double > values) |
Computes linear interpolation. | |
void | LagrangePolynomialInterpWeights (double x, std::span< const double > xs, std::vector< double > &weights) |
Computes weights for Lagrange polynomial interpolation. | |
void | LagrangePolynomialDerivWeights (double x, std::span< const double > xs, std::vector< double > &weights) |
Computes weights for Lagrange polynomial derivative interpolation. | |
double | LegendrePolynomial (const size_t n, const double x) |
Calculates the value of Legendre polynomial of order n at point x. | |
std::vector< double > | GenerateGLLPoints (const size_t order) |
Generates Gauss-Lobatto-Legendre (GLL) points for spectral element discretization. | |
std::vector< double > | MapGeometricLocations (std::span< const double > geom_locations) |
Maps input geometric locations -> normalized domain using linear mapping. | |
std::vector< std::vector< double > > | ComputeShapeFunctionValues (std::span< const double > input_points, std::span< const double > output_points) |
Computes shape function matrices ϕg relating points ξb to ξg At least two input points are required and it is assumed that there are more output points than input points. | |
std::vector< std::vector< double > > | ComputeShapeFunctionDerivatives (std::span< const double > input_points, std::span< const double > output_points) |
Computes shape function derivatives dϕg relating points ξb to ξg At least two input points are required and it is assumed that there are more output points than input points. | |
std::vector< std::array< double, 3 > > | PerformLeastSquaresFitting (size_t p, std::span< const std::vector< double > > shape_functions, std::span< const std::array< double, 3 > > points_to_fit) |
Performs least squares fitting to determine interpolation coefficients. | |
template<typename Matrix > | |
KOKKOS_INLINE_FUNCTION void | AX_Matrix (const Matrix &A, const Matrix &AX_A) |
Computes AX(A) of a square matrix. | |
template<typename Matrix , typename Vector > | |
KOKKOS_INLINE_FUNCTION void | AxialVectorOfMatrix (const Matrix &m, const Vector &v) |
Computes the axial vector (also known as the vector representation) of a 3x3 skew-symmetric matrix. | |
std::array< std::array< double, 6 >, 6 > | RotateMatrix6 (const std::array< std::array< double, 6 >, 6 > &m, const std::array< double, 4 > &q) |
std::vector< std::array< double, 3 > > | ProjectPointsToTargetPolynomial (size_t num_inputs, size_t num_outputs, std::span< const std::array< double, 3 > > input_points) |
Projects 3D points from a given (lower) polynomial representation to a target (higher) polynomial representation. | |
template<typename Quaternion , typename RotationMatrix > | |
KOKKOS_INLINE_FUNCTION void | QuaternionToRotationMatrix (const Quaternion &q, const RotationMatrix &R) |
Converts a 4x1 quaternion to a 3x3 rotation matrix and returns the result. | |
std::array< std::array< double, 3 >, 3 > | QuaternionToRotationMatrix (const std::array< double, 4 > &q) |
Converts a 4x1 quaternion to a 3x3 rotation matrix and returns the result. | |
std::array< double, 4 > | RotationMatrixToQuaternion (const std::array< std::array< double, 3 >, 3 > &m) |
Converts a 3x3 rotation matrix to a 4x1 quaternion and returns the result, see https://www.euclideanspace.com/maths/geometry/rotations/conversions/matrixToQuaternion/ for implementation details. | |
template<typename Quaternion , typename View1 , typename View2 > | |
KOKKOS_INLINE_FUNCTION void | RotateVectorByQuaternion (const Quaternion &q, const View1 &v, const View2 &v_rot) |
Rotates provided vector by provided unit quaternion and returns the result. | |
std::array< double, 3 > | RotateVectorByQuaternion (std::span< const double, 4 > q, std::span< const double, 3 > v) |
Rotates provided vector by provided unit quaternion and returns the result. | |
template<typename Quaternion , typename Matrix > | |
KOKKOS_INLINE_FUNCTION void | QuaternionDerivative (const Quaternion &q, const Matrix &m) |
Computes the derivative of a quaternion and stores the result in a 3x4 matrix. | |
template<typename QuaternionInput , typename QuaternionOutput > | |
KOKKOS_INLINE_FUNCTION void | QuaternionInverse (const QuaternionInput &q_in, const QuaternionOutput &q_out) |
Computes the inverse of a quaternion. | |
std::array< double, 4 > | QuaternionInverse (std::span< const double, 4 > quaternion) |
Computes the inverse of a quaternion. | |
template<typename Quaternion1 , typename Quaternion2 , typename QuaternionN > | |
KOKKOS_INLINE_FUNCTION void | QuaternionCompose (const Quaternion1 &q1, const Quaternion2 &q2, QuaternionN &qn) |
Composes (i.e. multiplies) two quaternions and stores the result in a third quaternion. | |
std::array< double, 4 > | QuaternionCompose (const std::array< double, 4 > &q1, const std::array< double, 4 > &q2) |
Composes (i.e. multiplies) two quaternions and returns the result. | |
template<typename Vector , typename Quaternion > | |
KOKKOS_INLINE_FUNCTION void | RotationVectorToQuaternion (const Vector &phi, const Quaternion &quaternion) |
Returns a 4-D quaternion from provided 3-D rotation vector, i.e. the exponential map. | |
template<typename Quaternion , typename Vector > | |
KOKKOS_INLINE_FUNCTION void | QuaternionToRotationVector (const Quaternion &quaternion, const Vector &phi) |
Returns a 3-D rotation vector from provided 4-D quaternion, i.e. the logarithmic map. | |
std::array< double, 3 > | QuaternionToRotationVector (const std::array< double, 4 > &quaternion) |
Returns a 3-D rotation vector from provided 4-D quaternion. | |
std::array< double, 4 > | RotationVectorToQuaternion (const std::array< double, 3 > &phi) |
Returns a 4-D quaternion from provided 3-D rotation vector, i.e. the exponential map. | |
KOKKOS_INLINE_FUNCTION Kokkos::Array< double, 4 > | NormalizeQuaternion (const Kokkos::Array< double, 4 > &q) |
Normalizes a quaternion to ensure it is a unit quaternion. | |
std::array< double, 4 > | TangentTwistToQuaternion (const std::array< double, 3 > &tangent, const double twist) |
Returns a 4-D quaternion from provided tangent vector and twist (degrees) about tangent. | |
std::array< double, 4 > | AxisAngleToQuaternion (const std::array< double, 3 > &axis, double angle) |
Returns a quaternion from the axis vector and angle (radians) | |
bool | IsIdentityQuaternion (const std::array< double, 4 > &q, double tolerance=1e-12) |
Checks if a quaternion is approximately the identity quaternion [1, 0, 0, 0]. | |
template<typename VectorType , typename MatrixType > | |
KOKKOS_INLINE_FUNCTION void | VecTilde (const VectorType &vector, const MatrixType &matrix) |
Converts a 3x1 vector to a 3x3 skew-symmetric matrix and returns the result. | |
template<typename AVectorType , typename BVectorType > | |
KOKKOS_INLINE_FUNCTION double | DotProduct (const AVectorType &a, const BVectorType &b) |
Calculate the dot product between two vector views. | |
constexpr double | DotProduct (const std::array< double, 3 > &a, const std::array< double, 3 > &b) |
Calculate the dot product between two vector views. | |
template<typename VectorType > | |
KOKKOS_INLINE_FUNCTION void | CrossProduct (const VectorType &a, const VectorType &b, const VectorType &c) |
Calculate the cross product between two vector views. | |
constexpr std::array< double, 3 > | CrossProduct (std::span< const double, 3 > a, std::span< const double, 3 > b) |
Calculate the cross product between two vectors. | |
constexpr double | Norm (const std::array< double, 3 > &v) |
Calculate the norm of a given vector. | |
constexpr std::array< double, 3 > | UnitVector (const std::array< double, 3 > &v) |
UnitVector returns the unit vector of the given vector. | |
Function Documentation
◆ AX_Matrix()
KOKKOS_INLINE_FUNCTION void kynema::math::AX_Matrix | ( | const Matrix & | A, |
const Matrix & | AX_A | ||
) |
Computes AX(A) of a square matrix.
AX(A) = tr(A)/2 * I - A/2, where I is the identity matrix
- Parameters
-
A Input square matrix AX_A Output matrix containing the result
◆ AxialVectorOfMatrix()
KOKKOS_INLINE_FUNCTION void kynema::math::AxialVectorOfMatrix | ( | const Matrix & | m, |
const Vector & | v | ||
) |
Computes the axial vector (also known as the vector representation) of a 3x3 skew-symmetric matrix.
The axial vector is defined as [w₁, w₂, w₃]ᵀ where: w₁ = (m₃₂ - m₂₃)/2 w₂ = (m₁₃ - m₃₁)/2 w₃ = (m₂₁ - m₁₂)/2
- Parameters
-
m Input 3x3 rotation matrix v Output vector to store the result
- Precondition
- Matrix m must be 3x3 and vector v must have size 3
◆ AxisAngleToQuaternion()
|
inline |
Returns a quaternion from the axis vector and angle (radians)
◆ ComputeShapeFunctionDerivatives()
|
inline |
Computes shape function derivatives dϕg relating points ξb to ξg At least two input points are required and it is assumed that there are more output points than input points.
- Parameters
-
input_points Input points, ξb, in [-1, 1] output_points Output points, ξg, in [-1, 1]
- Returns
- Shape function derivative matrix
◆ ComputeShapeFunctionValues()
|
inline |
Computes shape function matrices ϕg relating points ξb to ξg At least two input points are required and it is assumed that there are more output points than input points.
- Parameters
-
input_points Input points, ξb, in [-1, 1] output_points Output points, ξg, in [-1, 1]
- Returns
- Shape function matrix
◆ CrossProduct() [1/2]
KOKKOS_INLINE_FUNCTION void kynema::math::CrossProduct | ( | const VectorType & | a, |
const VectorType & | b, | ||
const VectorType & | c | ||
) |
Calculate the cross product between two vector views.
◆ CrossProduct() [2/2]
|
constexpr |
Calculate the cross product between two vectors.
◆ DotProduct() [1/2]
KOKKOS_INLINE_FUNCTION double kynema::math::DotProduct | ( | const AVectorType & | a, |
const BVectorType & | b | ||
) |
Calculate the dot product between two vector views.
◆ DotProduct() [2/2]
|
constexpr |
Calculate the dot product between two vector views.
◆ GenerateGLLPoints()
|
inline |
Generates Gauss-Lobatto-Legendre (GLL) points for spectral element discretization.
Computes the GLL points, i.e. roots of the Legendre polynomial, using Newton-Raphson iteration. GLL points are optimal interpolation nodes for spectral methods.
- Parameters
-
order Order of the polynomial interpolation (must be >= 1)
- Returns
- Vector of GLL points sorted in ascending order, size = order + 1
- Exceptions
-
std::invalid_argument if order < 1 std::runtime_error if Newton-Raphson iteration fails to converge
◆ GetGllLocations()
|
inline |
◆ GetGllWeights()
|
inline |
◆ IsIdentityQuaternion()
|
inline |
Checks if a quaternion is approximately the identity quaternion [1, 0, 0, 0].
- Parameters
-
q The quaternion to check tolerance The tolerance for the comparison (default: 1e-12)
- Returns
- true if the quaternion is approximately the identity quaternion, false otherwise
◆ LagrangePolynomialDerivWeights()
|
inline |
Computes weights for Lagrange polynomial derivative interpolation.
- Parameters
-
x Evaluation point xs Interpolation nodes (sorted) weights Output: weights for Lagrange polynomial derivative interpolation
◆ LagrangePolynomialInterpWeights()
|
inline |
Computes weights for Lagrange polynomial interpolation.
- Parameters
-
x Evaluation point xs Interpolation nodes (sorted) weights Output: weights for Lagrange polynomial interpolation
◆ LegendrePolynomial()
|
inline |
Calculates the value of Legendre polynomial of order n at point x.
Uses the recurrence relation for Legendre polynomials: P_n(x) = ((2n-1)xP_{n-1}(x) - (n-1)P_{n-2}(x))/n Reference: Deville et al. (2002) "High-Order Methods for Incompressible Fluid Flow" DOI: https://doi.org/10.1017/CBO9780511546792, Eq. B.1.15, p.446
- Parameters
-
n Order of the Legendre polynomial (n >= 0) x Point at which to evaluate the polynomial, typically in [-1,1]
- Returns
- Value of the nth order Legendre polynomial at x
◆ LinearInterp()
|
inline |
Computes linear interpolation.
- Parameters
-
x Evaluation point xs Value locations values Values at given locations
- Returns
- Interpolated value at evaluation point
◆ LinearInterpWeights()
|
inline |
Computes weights for linear interpolation.
- Parameters
-
x Evaluation point xs Interpolation nodes (sorted) weights Output: weights for linear interpolation
◆ MapGeometricLocations()
|
inline |
Maps input geometric locations -> normalized domain using linear mapping.
- Parameters
-
geom_locations Input geometric locations (typically in domain [0, 1]), sorted in ascending order
- Returns
- std::vector<double> Mapped/normalized evaluation points in domain [-1, 1]
◆ Norm()
|
constexpr |
Calculate the norm of a given vector.
◆ NormalizeQuaternion()
KOKKOS_INLINE_FUNCTION Kokkos::Array< double, 4 > kynema::math::NormalizeQuaternion | ( | const Kokkos::Array< double, 4 > & | q | ) |
Normalizes a quaternion to ensure it is a unit quaternion.
If the length of the quaternion is zero, it returns a default unit quaternion. Otherwise, it normalizes the quaternion and returns the result.
- Parameters
-
q The input quaternion as a Kokkos::Array<double, 4>
- Returns
- Kokkos::Array<double, 4> The normalized quaternion
◆ PerformLeastSquaresFitting()
|
inline |
Performs least squares fitting to determine interpolation coefficients.
Performs least squares fitting to determine interpolation coefficients by solving a dense linear system [A][X] = [B], where [A] is the shape function matrix (p x n), [B] is the input points (n x 3), and [X] is the interpolation coefficients (p x 3)
- Parameters
-
p Number of points representing the polynomial of order p-1 shape_functions Shape function matrix (p x n) points_to_fit x,y,z coordinates of the points to fit (n x 3)
- Returns
- Interpolation coefficients (p x 3)
◆ ProjectPointsToTargetPolynomial()
|
inline |
Projects 3D points from a given (lower) polynomial representation to a target (higher) polynomial representation.
This function maps a set of 3D points defined at nodes of a polynomial of order source_order to corresponding points at nodes of a polynomial of order target_order (typically higher than the source order) using Least-Squares Finite Element (LSFE) shape functions.
Primary use case: The primary application of this function is to increase the number of points from a lower-order geometric representation to the higher-order representation required for spectral finite element analysis. This enables the use of high-order methods while allowing geometry to be defined with fewer points initially.
- Parameters
-
num_inputs Number of points in the source polynomial representation num_outputs Number of points in the target polynomial representation input_points 3D coordinates of points in the source representation
- Returns
- std::vector<std::array<double, 3>>
- Coordinates of the projected 3D points at the target polynomial nodes
◆ QuaternionCompose() [1/2]
KOKKOS_INLINE_FUNCTION void kynema::math::QuaternionCompose | ( | const Quaternion1 & | q1, |
const Quaternion2 & | q2, | ||
QuaternionN & | qn | ||
) |
Composes (i.e. multiplies) two quaternions and stores the result in a third quaternion.
◆ QuaternionCompose() [2/2]
|
inline |
Composes (i.e. multiplies) two quaternions and returns the result.
◆ QuaternionDerivative()
KOKKOS_INLINE_FUNCTION void kynema::math::QuaternionDerivative | ( | const Quaternion & | q, |
const Matrix & | m | ||
) |
Computes the derivative of a quaternion and stores the result in a 3x4 matrix.
◆ QuaternionInverse() [1/2]
KOKKOS_INLINE_FUNCTION void kynema::math::QuaternionInverse | ( | const QuaternionInput & | q_in, |
const QuaternionOutput & | q_out | ||
) |
Computes the inverse of a quaternion.
◆ QuaternionInverse() [2/2]
|
inline |
Computes the inverse of a quaternion.
◆ QuaternionToRotationMatrix() [1/2]
KOKKOS_INLINE_FUNCTION void kynema::math::QuaternionToRotationMatrix | ( | const Quaternion & | q, |
const RotationMatrix & | R | ||
) |
Converts a 4x1 quaternion to a 3x3 rotation matrix and returns the result.
◆ QuaternionToRotationMatrix() [2/2]
|
inline |
Converts a 4x1 quaternion to a 3x3 rotation matrix and returns the result.
◆ QuaternionToRotationVector() [1/2]
KOKKOS_INLINE_FUNCTION void kynema::math::QuaternionToRotationVector | ( | const Quaternion & | quaternion, |
const Vector & | phi | ||
) |
Returns a 3-D rotation vector from provided 4-D quaternion, i.e. the logarithmic map.
◆ QuaternionToRotationVector() [2/2]
|
inline |
Returns a 3-D rotation vector from provided 4-D quaternion.
◆ RotateMatrix6()
|
inline |
◆ RotateVectorByQuaternion() [1/2]
KOKKOS_INLINE_FUNCTION void kynema::math::RotateVectorByQuaternion | ( | const Quaternion & | q, |
const View1 & | v, | ||
const View2 & | v_rot | ||
) |
Rotates provided vector by provided unit quaternion and returns the result.
◆ RotateVectorByQuaternion() [2/2]
|
inline |
Rotates provided vector by provided unit quaternion and returns the result.
◆ RotationMatrixToQuaternion()
|
inline |
Converts a 3x3 rotation matrix to a 4x1 quaternion and returns the result, see https://www.euclideanspace.com/maths/geometry/rotations/conversions/matrixToQuaternion/ for implementation details.
◆ RotationVectorToQuaternion() [1/2]
|
inline |
Returns a 4-D quaternion from provided 3-D rotation vector, i.e. the exponential map.
◆ RotationVectorToQuaternion() [2/2]
KOKKOS_INLINE_FUNCTION void kynema::math::RotationVectorToQuaternion | ( | const Vector & | phi, |
const Quaternion & | quaternion | ||
) |
Returns a 4-D quaternion from provided 3-D rotation vector, i.e. the exponential map.
◆ TangentTwistToQuaternion()
|
inline |
Returns a 4-D quaternion from provided tangent vector and twist (degrees) about tangent.
◆ UnitVector()
|
constexpr |
UnitVector returns the unit vector of the given vector.
◆ VecTilde()
KOKKOS_INLINE_FUNCTION void kynema::math::VecTilde | ( | const VectorType & | vector, |
const MatrixType & | matrix | ||
) |
Converts a 3x1 vector to a 3x3 skew-symmetric matrix and returns the result.
Generated by