Description

Collaboration diagram for Linear algebra:

Classes

class  chrono::ChMatrix33< Real >
 Definition of a 3x3 fixed size matrix to represent 3D rotations and inertia tensors. More...
 
class  chrono::ChFpMatrix34< Real >
 Special MBD 3x4 matrix [Fp(q)], as in [Fp(q)] * [Fm(q)]' = [A(q)]. More...
 
class  chrono::ChFmMatrix34< Real >
 Special MBD 3x4 matrix [Fm(q)], as in [Fp(q)] * [Fm(q)]' = [A(q)]. More...
 
class  chrono::ChGlMatrix34< Real >
 Special MBD 3x4 matrix [Gl(q)], as in local angular speed conversion. More...
 
class  chrono::ChGwMatrix34< Real >
 Special MBD 3x4 matrix [Gw(q)], as in absolute angular speed conversion. More...
 
class  chrono::ChStarMatrix33< Real >
 Special MBD 3x3 "star" matrix, , representing vector cross products. More...
 
class  chrono::ChStarMatrix44< Real >
 Special MBD 4x4 "star" matrix, representing quaternion cross product. More...
 
class  chrono::ChSparsityPatternLearner
 Utility class for extracting sparsity patter from a sparse matrix. More...
 

Typedefs

template<typename T = double>
using chrono::ChMatrixDynamic = Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >
 Dense matrix with dynamic size (i.e., with unknown at compile time) and row-major storage. More...
 
template<typename T , int M, int N>
using chrono::ChMatrixNM = Eigen::Matrix< T, M, N, Eigen::RowMajor >
 Dense matrix with fixed size (known at compile time) and row-major storage. More...
 
template<typename T = double>
using chrono::ChMatrixDynamic_col = Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >
 Dense matrix with dynamic size (i.e., with unknown at compile time) and column-major storage. More...
 
template<typename T , int M, int N>
using chrono::ChMatrixNM_col = Eigen::Matrix< T, M, N, Eigen::ColMajor >
 Dense matrix with fixed size (known at compile time) and column-major storage. More...
 
template<typename T = double>
using chrono::ChVectorDynamic = Eigen::Matrix< T, Eigen::Dynamic, 1, Eigen::ColMajor >
 Column vector with dynamic size (i.e., with size unknown at compile time). More...
 
template<typename T = double>
using chrono::ChRowVectorDynamic = Eigen::Matrix< T, 1, Eigen::Dynamic, Eigen::RowMajor >
 Row vector with dynamic size (i.e., with size unknown at compile time). More...
 
template<typename T , int N>
using chrono::ChVectorN = Eigen::Matrix< T, N, 1 >
 Column vector with fixed size (known at compile time). More...
 
template<typename T , int N>
using chrono::ChRowVectorN = Eigen::Matrix< T, 1, N, Eigen::RowMajor >
 Row vector with fixed size (known at compile time). More...
 
template<typename T = double>
using chrono::ChArray = Eigen::Array< T, Eigen::Dynamic, 1, Eigen::ColMajor >
 General-purpose column array with dynamic size. More...
 
template<typename T , int N>
using chrono::ChArrayN = Eigen::Array< T, N, 1 >
 Column array with fixed size (known at compile time). More...
 
using chrono::ChMatrixRef = Eigen::Ref< Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > >
 Reference to a dense matrix expression, with double coefficients. More...
 
using chrono::ChMatrixConstRef = const Eigen::Ref< const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > > &
 Constant reference to a dense matrix expression, with double coefficients. More...
 
using chrono::ChVectorRef = Eigen::Ref< Eigen::Matrix< double, Eigen::Dynamic, 1, Eigen::ColMajor > >
 Reference to a column vector expression, with double coefficients. More...
 
using chrono::ChVectorConstRef = const Eigen::Ref< const Eigen::Matrix< double, Eigen::Dynamic, 1, Eigen::ColMajor > > &
 Constant reference to a column vector expression, with double coefficients. More...
 
using chrono::ChRowVectorRef = Eigen::Ref< Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor > >
 Reference to a row vector expression, with double coefficients. More...
 
using chrono::ChRowVectorConstRef = const Eigen::Ref< const Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor > > &
 Constant reference to a row vector expression, with double coefficients. More...
 
template<typename T = double>
using chrono::ChArrayRef = Eigen::Ref< Eigen::Array< T, Eigen::Dynamic, Eigen::ColMajor > > &
 Reference to an array expression, templated by coefficient type. More...
 
template<typename T = double>
using chrono::ChArrayConstRef = const Eigen::Ref< const Eigen::Array< T, Eigen::Dynamic, 1, Eigen::ColMajor > > &
 Constant reference to an array expression, templated by coefficient type. More...
 
using chrono::ChSparseMatrix = Eigen::SparseMatrix< double, Eigen::RowMajor, int >
 Sparse matrix representation. More...
 
template<typename T = double>
using chrono::ChMatrix34 = Eigen::Matrix< T, 3, 4, Eigen::RowMajor >
 3x4 dense matrix (fixed-size, row-major ordering)
 
template<typename T = double>
using chrono::ChMatrix43 = Eigen::Matrix< T, 4, 3, Eigen::ColMajor >
 4x3 dense matrix (fixed-size, row-major ordering)
 
template<typename T = double>
using chrono::ChMatrix44 = Eigen::Matrix< T, 4, 4, Eigen::RowMajor >
 4x4 dense matrix (fixed-size, row-major ordering)
 

Functions

void chrono::StreamOutDenseMatlabFormat (ChMatrixConstRef A, ChStreamOutAscii &stream)
 Serialization of a dense matrix or vector into an ASCII stream (e.g. a file) in Matlab format.
 
void chrono::StreamInDenseMatlabFormat (const std::string &filename, ChMatrixDynamic<> &matr, char delim=',')
 Parse numeric data from a file (eg. csv) and store into dense matrix.
 
void chrono::PasteMatrix (ChSparseMatrix &matrTo, ChMatrixConstRef matrFrom, int insrow, int inscol, bool overwrite=true)
 Paste a given matrix into a sparse matrix at position (insrow, inscol). More...
 
void chrono::StreamOutSparseMatlabFormat (ChSparseMatrix &matr, ChStreamOutAscii &mstream)
 Serialization of a sparse matrix to an ASCI stream (e.g., a file) in Matlab sparse matrix format. More...
 
void chrono::StreamOut (ChSparseMatrix &matr, ChStreamOutAscii &stream)
 Serialization of a sparse matrix to an ASCII stream (for debugging; only the top-left 8x8 corner is printed).
 
template<typename T = double>
ChVectorDynamic< T > chrono::SliceVector (ChVectorConstRef v, ChArrayConstRef< int > indices)
 Utility function for slicing a vector based on an array of indices. More...
 
template<typename Real >
ChVector< Real > chrono::operator* (const Eigen::Transpose< Eigen::Matrix< Real, 3, 3, Eigen::RowMajor >> &A, const ChVector< Real > &v)
 Multiply a transposed 3x3 matrix with a vector.
 
template<typename Real >
ChVector< Real > chrono::operator* (const Eigen::Transpose< const Eigen::Matrix< Real, 3, 3, Eigen::RowMajor >> &A, const ChVector< Real > &v)
 Multiply a transposed const 3x3 matrix with a vector.
 
template<class Real >
ChMatrix33< Real > chrono::TensorProduct (const ChVector< Real > &vA, const ChVector< Real > &vB)
 Return the outer product (a 3x3 matrix) of two vectors.
 
template<typename T , typename U >
ChVector< T > chrono::operator* (const ChMatrix34< T > &A, const ChQuaternion< U > &q)
 Multiply a 3x4 matrix with a quaternion and return a 3d vector.
 
template<typename T , typename U >
ChQuaternion< T > chrono::operator* (const ChMatrix43< T > &A, const ChVector< U > &v)
 Multiply a 4x3 matrix with a 3d vector and return a quaternion.
 
template<typename T , typename U >
ChQuaternion< T > chrono::operator* (const Eigen::Transpose< Eigen::Matrix< T, 3, 4, Eigen::RowMajor >> &A, const ChVector< U > &v)
 Multiply the transpose of a 3x4 matrix with a 3d vector and return a quaternion.
 
template<typename T , typename U >
ChQuaternion< T > chrono::operator* (const ChMatrix44< T > &A, const ChQuaternion< U > &q)
 Multiply a 4x4 matrix with a quaternion and return a quaternion.
 
 chrono::ChMatrix33< Real >::ChMatrix33 (const ChQuaternion< Real > &q)
 Construct a 3x3 rotation matrix from the given quaternion.
 
 chrono::ChMatrix33< Real >::ChMatrix33 (Real val)
 Construct a diagonal matrix with the specified value on the diagonal.
 
 chrono::ChMatrix33< Real >::ChMatrix33 (const ChVector< Real > &v)
 Construct a diagonal matrix with the specified values on the diagonal.
 
 chrono::ChMatrix33< Real >::ChMatrix33 (const ChVector<> &diag, const ChVector<> &off_diag)
 Construct a symmetric 3x3 matrix with the specified vectors for the diagonal and off-digonal elements. More...
 
 chrono::ChMatrix33< Real >::ChMatrix33 (Real angle, const ChVector<> &axis)
 Construct a 3x3 rotation matrix from an angle and a rotation axis. More...
 
 chrono::ChMatrix33< Real >::ChMatrix33 (const ChVector<> &X, const ChVector<> &Y, const ChVector<> &Z)
 Construct a 3x3 matrix with the given vectors as columns. More...
 
ChVector< Real > chrono::ChMatrix33< Real >::operator* (const ChVector< Real > &v) const
 Multiply this matrix by a 3d vector.
 
void chrono::ChMatrix33< Real >::Set_A_quaternion (const ChQuaternion< Real > &quat)
 Fill this 3x3 matrix as a rotation matrix, given a unit quaternion.
 
void chrono::ChMatrix33< Real >::Set_A_Eulero (const ChVector< Real > &angles)
 Fill this 3x3 matrix as a rotation matrix, given three Euler angles.
 
void chrono::ChMatrix33< Real >::Set_A_Cardano (const ChVector< Real > &angles)
 Fill this 3x3 matrix as a rotation matrix, given three Cardano angles.
 
void chrono::ChMatrix33< Real >::Set_A_Hpb (const ChVector< Real > &angles)
 Fill this 3x3 matrix as a rotation matrix, given three head, pitch, banking angles.
 
void chrono::ChMatrix33< Real >::Set_A_Rxyz (const ChVector< Real > &xyz)
 Fill this 3x3 matrix as a rotation matrix, given three angles of consecutive rotations about x,y,z axis.
 
void chrono::ChMatrix33< Real >::Set_A_Rodriguez (const ChVector< Real > &rod)
 Fill this 3x3 matrix as a rotation matrix, given three Rodriguez parameters.
 
void chrono::ChMatrix33< Real >::Set_A_axis (const ChVector< Real > &X, const ChVector< Real > &Y, const ChVector< Real > &Z)
 Fill this 3x3 matrix as a rotation matrix, given the three versors X,Y,Z of the basis.
 
void chrono::ChMatrix33< Real >::Set_A_Xdir (const ChVector< Real > &Xdir, const ChVector< Real > &Vsingular=ChVector< Real >(0, 1, 0))
 Fill this 3x3 matrix as a rotation matrix with the X axis along the provided direction. More...
 
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Eulero () const
 Return the Eulero angles. More...
 
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Cardano () const
 Return the Cardano angles. More...
 
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Hpb () const
 Return the head-pitch-banking angles. More...
 
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Rxyz () const
 Return the angles for consecutive rotations on x,y,z axes. More...
 
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Rodriguez () const
 Return the Rodriguez parameters. More...
 
ChQuaternion< Real > chrono::ChMatrix33< Real >::Get_A_quaternion () const
 Return the corresponding unit quaternion. More...
 
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Xaxis () const
 Return the versor of X axis.
 
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Yaxis () const
 Return the versor of Y axis.
 
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Zaxis () const
 Return the versor of Z axis.
 
ChVector< Real > chrono::ChMatrix33< Real >::GetAx () const
 Assuming this matrix is a rotation matrix, get Ax vector.
 
void chrono::ChMatrix33< Real >::SelfAdjointEigenSolve (ChMatrix33< Real > &evec, ChVectorN< Real, 3 > &evals) const
 Compute eigenvectors and eigenvalues. More...
 

Typedef Documentation

◆ ChArray

template<typename T = double>
using chrono::ChArray = typedef Eigen::Array<T, Eigen::Dynamic, 1, Eigen::ColMajor>

General-purpose column array with dynamic size.

This class provides easy-access to coefficient-wise operations.

◆ ChArrayConstRef

template<typename T = double>
using chrono::ChArrayConstRef = typedef const Eigen::Ref<const Eigen::Array<T, Eigen::Dynamic, 1, Eigen::ColMajor> >&

Constant reference to an array expression, templated by coefficient type.

This allows writing non-template functions that can accept either a ChArray or a ChArrayN.

◆ ChArrayN

template<typename T , int N>
using chrono::ChArrayN = typedef Eigen::Array<T, N, 1>

Column array with fixed size (known at compile time).

A ChArrayN is templated by the type of its coefficients and the number of elements.

◆ ChArrayRef

template<typename T = double>
using chrono::ChArrayRef = typedef Eigen::Ref<Eigen::Array<T, Eigen::Dynamic, Eigen::ColMajor> >&

Reference to an array expression, templated by coefficient type.

This allows writing non-template functions that can accept either a ChArrayDynamic or a ChArrayN.

◆ ChMatrixConstRef

using chrono::ChMatrixConstRef = typedef const Eigen::Ref<const Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> >&

Constant reference to a dense matrix expression, with double coefficients.

This allows writing non-template functions that can accept either a ChMatrixDynamic or a ChMatrixNM.

◆ ChMatrixDynamic

template<typename T = double>
using chrono::ChMatrixDynamic = typedef Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>

Dense matrix with dynamic size (i.e., with unknown at compile time) and row-major storage.

A ChMatrixDynamic is templated by the coefficient type (default: double).

◆ ChMatrixDynamic_col

template<typename T = double>
using chrono::ChMatrixDynamic_col = typedef Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor>

Dense matrix with dynamic size (i.e., with unknown at compile time) and column-major storage.

A ChMatrixDynamic_col is templated by the coefficient type (default: double).

◆ ChMatrixNM

template<typename T , int M, int N>
using chrono::ChMatrixNM = typedef Eigen::Matrix<T, M, N, Eigen::RowMajor>

Dense matrix with fixed size (known at compile time) and row-major storage.

A ChMatrixNM is templated by the coefficient type and by the matrix dimensions (number of rows and columns).

◆ ChMatrixNM_col

template<typename T , int M, int N>
using chrono::ChMatrixNM_col = typedef Eigen::Matrix<T, M, N, Eigen::ColMajor>

Dense matrix with fixed size (known at compile time) and column-major storage.

A ChMatrixNM_col is templated by the coefficient type and by the matrix dimensions (number of rows and columns).

◆ ChMatrixRef

using chrono::ChMatrixRef = typedef Eigen::Ref<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> >

Reference to a dense matrix expression, with double coefficients.

This allows writing non-template functions that can accept either a ChMatrixDynamic or a ChMatrixNM.

◆ ChRowVectorConstRef

using chrono::ChRowVectorConstRef = typedef const Eigen::Ref<const Eigen::Matrix<double, 1, Eigen::Dynamic, Eigen::RowMajor> >&

Constant reference to a row vector expression, with double coefficients.

This allows writing non-template functions that can accept either a ChRowVectorDynamic or a CVectorN.

◆ ChRowVectorDynamic

template<typename T = double>
using chrono::ChRowVectorDynamic = typedef Eigen::Matrix<T, 1, Eigen::Dynamic, Eigen::RowMajor>

Row vector with dynamic size (i.e., with size unknown at compile time).

A ChRowVectorDynamic is templated by the type of its coefficients (default: double).

◆ ChRowVectorN

template<typename T , int N>
using chrono::ChRowVectorN = typedef Eigen::Matrix<T, 1, N, Eigen::RowMajor>

Row vector with fixed size (known at compile time).

A ChRowVectorN is templated by the type of its coefficients and the number of elements.

◆ ChRowVectorRef

using chrono::ChRowVectorRef = typedef Eigen::Ref<Eigen::Matrix<double, 1, Eigen::Dynamic, Eigen::RowMajor> >

Reference to a row vector expression, with double coefficients.

This allows writing non-template functions that can accept either a ChRowVectorDynamic or a CVectorN.

◆ ChSparseMatrix

using chrono::ChSparseMatrix = typedef Eigen::SparseMatrix<double, Eigen::RowMajor, int>

Sparse matrix representation.

A ChSparseMatrix is an Eigen SparseMatrix with double coefficients, row-major storage order, and int indices.

◆ ChVectorConstRef

using chrono::ChVectorConstRef = typedef const Eigen::Ref<const Eigen::Matrix<double, Eigen::Dynamic, 1, Eigen::ColMajor> >&

Constant reference to a column vector expression, with double coefficients.

This allows writing non-template functions that can accept either a ChVectorDynamic or a ChRowVectorN.

◆ ChVectorDynamic

template<typename T = double>
using chrono::ChVectorDynamic = typedef Eigen::Matrix<T, Eigen::Dynamic, 1, Eigen::ColMajor>

Column vector with dynamic size (i.e., with size unknown at compile time).

A ChVectorDynamic is templated by the type of its coefficients (default: double).

◆ ChVectorN

template<typename T , int N>
using chrono::ChVectorN = typedef Eigen::Matrix<T, N, 1>

Column vector with fixed size (known at compile time).

A ChVectorN is templated by the type of its coefficients and the number of elements.

◆ ChVectorRef

using chrono::ChVectorRef = typedef Eigen::Ref<Eigen::Matrix<double, Eigen::Dynamic, 1, Eigen::ColMajor> >

Reference to a column vector expression, with double coefficients.

This allows writing non-template functions that can accept either a ChVectorDynamic or a ChVectorN.

Function Documentation

◆ ChMatrix33() [1/3]

template<typename Real>
chrono::ChMatrix33< Real >::ChMatrix33 ( const ChVector<> &  diag,
const ChVector<> &  off_diag 
)

Construct a symmetric 3x3 matrix with the specified vectors for the diagonal and off-digonal elements.

The off-diagonal vector is assumed to contain the elements A(0,1), A(0,2), A(1,2) in this order.

◆ ChMatrix33() [2/3]

template<typename Real>
chrono::ChMatrix33< Real >::ChMatrix33 ( const ChVector<> &  X,
const ChVector<> &  Y,
const ChVector<> &  Z 
)

Construct a 3x3 matrix with the given vectors as columns.

If the three vectors are mutually orthogonal unit vectors, the resulting matrix is a rotation matrix.

◆ ChMatrix33() [3/3]

template<typename Real>
chrono::ChMatrix33< Real >::ChMatrix33 ( Real  angle,
const ChVector<> &  axis 
)

Construct a 3x3 rotation matrix from an angle and a rotation axis.

Note that the axis direction must be normalized.

◆ Get_A_Cardano()

template<typename Real >
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Cardano ( ) const
inline

Return the Cardano angles.

Assumes that this is a rotation matrix.

◆ Get_A_Eulero()

template<typename Real >
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Eulero ( ) const
inline

Return the Eulero angles.

Assumes that this is a rotation matrix.

◆ Get_A_Hpb()

template<typename Real >
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Hpb ( ) const
inline

Return the head-pitch-banking angles.

Assumes that this is a rotation matrix.

◆ Get_A_quaternion()

template<typename Real >
ChQuaternion< Real > chrono::ChMatrix33< Real >::Get_A_quaternion ( ) const
inline

Return the corresponding unit quaternion.

Assumes that this is a rotation matrix.

◆ Get_A_Rodriguez()

template<typename Real >
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Rodriguez ( ) const
inline

Return the Rodriguez parameters.

Assumes that this is a rotation matrix.

◆ Get_A_Rxyz()

template<typename Real >
ChVector< Real > chrono::ChMatrix33< Real >::Get_A_Rxyz ( ) const
inline

Return the angles for consecutive rotations on x,y,z axes.

Assumes that this is a rotation matrix.

◆ PasteMatrix()

void chrono::PasteMatrix ( ChSparseMatrix matrTo,
ChMatrixConstRef  matrFrom,
int  insrow,
int  inscol,
bool  overwrite = true 
)
inline

Paste a given matrix into a sparse matrix at position (insrow, inscol).

The matrix matrFrom will be copied into matrTo[insrow : insrow + matrFrom.GetRows()][inscol : inscol + matrFrom.GetColumns()]

Parameters
[out]matrToThe output sparse matrix
[in]matrFromThe source matrix that will be copied
[in]insrowThe row index where the first element will be copied
[in]inscolThe column index where the first element will be copied
[in]overwriteIndicate if the copied elements will overwrite existing elements or be summed to them

◆ SelfAdjointEigenSolve()

template<typename Real>
void chrono::ChMatrix33< Real >::SelfAdjointEigenSolve ( ChMatrix33< Real > &  evec,
ChVectorN< Real, 3 > &  evals 
) const
inline

Compute eigenvectors and eigenvalues.

Note: only for self-adjoint matrices (e.g. inertia tensors).

◆ Set_A_Xdir()

template<typename Real>
void chrono::ChMatrix33< Real >::Set_A_Xdir ( const ChVector< Real > &  Xdir,
const ChVector< Real > &  Vsingular = ChVector<Real>(0, 1, 0) 
)
inline

Fill this 3x3 matrix as a rotation matrix with the X axis along the provided direction.

Uses the Gram-Schmidt orthonormalization. The optional argument Vsingular, together with Xdir, suggests the XY plane (as long as Xdir is not too close to lying in that plane, in which case a different direction is selected).

Parameters
XdirX axis
Vsingularsuggested Y axis

◆ SliceVector()

template<typename T = double>
ChVectorDynamic<T> chrono::SliceVector ( ChVectorConstRef  v,
ChArrayConstRef< int >  indices 
)

Utility function for slicing a vector based on an array of indices.

Return a new vector which only contains the elements with specified indices.

◆ StreamOutSparseMatlabFormat()

void chrono::StreamOutSparseMatlabFormat ( ChSparseMatrix matr,
ChStreamOutAscii mstream 
)
inline

Serialization of a sparse matrix to an ASCI stream (e.g., a file) in Matlab sparse matrix format.

Note that row and column indices start at 1.