Description
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 > | |
using | chrono::ChMatrix66 = ChMatrixNM< T, 6, 6 > |
Alias for a 6x6 matrix templated by coefficient type (row-major storage). | |
using | chrono::ChMatrix66d = ChMatrix66< double > |
Alias for a 6x6 matrix of doubles. | |
using | chrono::ChMatrix66f = ChMatrix66< float > |
Alias for a 6x6 matrix of floats. | |
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... | |
using | chrono::ChMatrix33d = ChMatrix33< double > |
Alias for a 3x3 matrix of doubles. | |
using | chrono::ChMatrix33f = ChMatrix33< float > |
Alias for a 3x3 matrix of floats. | |
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) | |
using | chrono::ChStarMatrix33d = ChStarMatrix33< double > |
Alias for a 3x3 star matrix of doubles. | |
Functions | |
void | chrono::PasteMatrix (ChSparseMatrix &matrTo, ChMatrixConstRef matrFrom, int start_row, int start_col, bool overwrite=true) |
Paste a given matrix into a sparse matrix at position (start_row, start_col). More... | |
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... | |
void | chrono::StreamOut (ChMatrixConstRef A, std::ostream &stream) |
Serialization of a dense matrix or vector into an ASCII stream (e.g. a file). | |
void | chrono::StreamOut (ChSparseMatrix &mat, std::ostream &stream, bool one_indexed=false) |
Serialization of a sparse matrix to an ASCII stream (e.g., a file) in COO sparse matrix format. More... | |
template<typename Real > | |
ChVector3< Real > | chrono::operator* (const Eigen::Transpose< Eigen::Matrix< Real, 3, 3, Eigen::RowMajor >> &A, const ChVector3< Real > &v) |
Multiply a transposed 3x3 matrix with a vector. | |
template<typename Real > | |
ChVector3< Real > | chrono::operator* (const Eigen::Transpose< const Eigen::Matrix< Real, 3, 3, Eigen::RowMajor >> &A, const ChVector3< Real > &v) |
Multiply a transposed const 3x3 matrix with a vector. | |
template<class Real > | |
ChMatrix33< Real > | chrono::TensorProduct (const ChVector3< Real > &vA, const ChVector3< Real > &vB) |
Return the outer product (a 3x3 matrix) of two vectors. | |
template<typename T , typename U > | |
ChVector3< 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 ChVector3< 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 ChVector3< 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 ChVector3< Real > &v) | |
Construct a diagonal matrix with the specified values on the diagonal. | |
chrono::ChMatrix33< Real >::ChMatrix33 (const ChVector3d &diag, const ChVector3d &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 ChVector3d &axis) | |
Construct a 3x3 rotation matrix from an angle and a rotation axis. More... | |
chrono::ChMatrix33< Real >::ChMatrix33 (const ChVector3d &X, const ChVector3d &Y, const ChVector3d &Z) | |
Construct a 3x3 matrix with the given vectors as columns. More... | |
ChVector3< Real > | chrono::ChMatrix33< Real >::operator* (const ChVector3< Real > &v) const |
Multiply this matrix by a 3d vector. | |
void | chrono::ChMatrix33< Real >::SetFromQuaternion (const ChQuaternion< Real > &q) |
Fill this 3x3 matrix as a rotation matrix, given a unit quaternion. | |
void | chrono::ChMatrix33< Real >::SetFromEulerAnglesZXZ (const ChVector3< Real > &angles) |
Fill this 3x3 matrix as a rotation matrix, given three Euler angles. | |
void | chrono::ChMatrix33< Real >::SetFromCardanAnglesZXY (const ChVector3< Real > &angles) |
Fill this 3x3 matrix as a rotation matrix, given three Cardano angles. | |
void | chrono::ChMatrix33< Real >::SetFromCardanAnglesZYX (const ChVector3< Real > &angles) |
Fill this 3x3 matrix as a rotation matrix, given three head, pitch, banking angles. | |
void | chrono::ChMatrix33< Real >::SetFromCardanAnglesXYZ (const ChVector3< Real > &angles) |
Fill this 3x3 matrix as a rotation matrix, given three angles of consecutive rotations about x,y,z axis. | |
void | chrono::ChMatrix33< Real >::SetFromRodriguesParameters (const ChVector3< Real > &r) |
Fill this 3x3 matrix as a rotation matrix, given three Rodrigues parameters. | |
void | chrono::ChMatrix33< Real >::SetFromDirectionAxes (const ChVector3< Real > &X, const ChVector3< Real > &Y, const ChVector3< Real > &Z) |
Fill this 3x3 matrix as a rotation matrix, given the three versors X,Y,Z of the basis. | |
void | chrono::ChMatrix33< Real >::SetFromAxisX (const ChVector3< Real > &x_dir, const ChVector3< Real > &y_sugg=ChVector3< Real >(0, 1, 0)) |
Fill this 3x3 matrix as a rotation matrix with the X axis along the provided direction. More... | |
void | chrono::ChMatrix33< Real >::SetFromAxisY (const ChVector3< Real > &y_dir, const ChVector3< Real > &z_sugg=ChVector3< Real >(0, 0, 1)) |
Fill this 3x3 matrix as a rotation matrix with the Y axis along the provided direction. More... | |
void | chrono::ChMatrix33< Real >::SetFromAxisZ (const ChVector3< Real > &z_dir, const ChVector3< Real > &x_sugg=ChVector3< Real >(1, 0, 0)) |
Fill this 3x3 matrix as a rotation matrix with the Z axis along the provided direction. More... | |
ChVector3< Real > | chrono::ChMatrix33< Real >::GetEulerAnglesZXZ () const |
Return the Euler angles. More... | |
ChVector3< Real > | chrono::ChMatrix33< Real >::GetCardanAnglesZXY () const |
Return the Cardano angles. More... | |
ChVector3< Real > | chrono::ChMatrix33< Real >::GetCardanAnglesZYX () const |
Return the head-pitch-banking angles. More... | |
ChVector3< Real > | chrono::ChMatrix33< Real >::GetCardanAnglesXYZ () const |
Return the angles for consecutive rotations on x,y,z axes. More... | |
ChVector3< Real > | chrono::ChMatrix33< Real >::GetRodriguesParameters () const |
Return the Rodrigues parameters. More... | |
ChQuaternion< Real > | chrono::ChMatrix33< Real >::GetQuaternion () const |
Return the corresponding unit quaternion. More... | |
ChVector3< Real > | chrono::ChMatrix33< Real >::GetAxisX () const |
Return the unit vector along the X axis. | |
ChVector3< Real > | chrono::ChMatrix33< Real >::GetAxisY () const |
Return the unit vector along the Y axis. | |
ChVector3< Real > | chrono::ChMatrix33< Real >::GetAxisZ () const |
Return the unit vector along the Z axis. | |
void | chrono::ChMatrix33< Real >::SelfAdjointEigenSolve (ChMatrix33< Real > &evec, ChVectorN< Real, 3 > &evals) const |
Compute eigenvectors and eigenvalues. More... | |
Typedef Documentation
◆ ChArray
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
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
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
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
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
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
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
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
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
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
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
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]
chrono::ChMatrix33< Real >::ChMatrix33 | ( | const ChVector3d & | diag, |
const ChVector3d & | 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]
chrono::ChMatrix33< Real >::ChMatrix33 | ( | const ChVector3d & | X, |
const ChVector3d & | Y, | ||
const ChVector3d & | 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]
chrono::ChMatrix33< Real >::ChMatrix33 | ( | Real | angle, |
const ChVector3d & | axis | ||
) |
Construct a 3x3 rotation matrix from an angle and a rotation axis.
Note that the axis direction must be normalized.
◆ GetCardanAnglesXYZ()
|
inline |
Return the angles for consecutive rotations on x,y,z axes.
Assumes that this is a rotation matrix.
◆ GetCardanAnglesZXY()
|
inline |
Return the Cardano angles.
Assumes that this is a rotation matrix.
◆ GetCardanAnglesZYX()
|
inline |
Return the head-pitch-banking angles.
Assumes that this is a rotation matrix.
◆ GetEulerAnglesZXZ()
|
inline |
Return the Euler angles.
Assumes that this is a rotation matrix.
◆ GetQuaternion()
|
inline |
Return the corresponding unit quaternion.
Assumes that this is a rotation matrix.
◆ GetRodriguesParameters()
|
inline |
Return the Rodrigues parameters.
Assumes that this is a rotation matrix.
◆ PasteMatrix()
|
inline |
Paste a given matrix into a sparse matrix at position (start_row, start_col).
The matrix matrFrom will be copied into matrTo[start_row : start_row + matrFrom.GetRows()][start_col : start_col + matrFrom.GetColumns()]
- Parameters
-
[out] matrTo The output sparse matrix [in] matrFrom The source matrix that will be copied [in] start_row The row index where the first element will be copied [in] start_col The column index where the first element will be copied [in] overwrite Indicate if the copied elements will overwrite existing elements or be summed to them
◆ SelfAdjointEigenSolve()
|
inline |
Compute eigenvectors and eigenvalues.
Note: only for self-adjoint matrices (e.g. inertia tensors).
◆ SetFromAxisX()
|
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 y_sugg, together with x_dir, suggests the XY plane (as long as y_sugg is not too close y_sugg, in which case a different direction is selected).
- Parameters
-
x_dir X axis y_sugg suggested Y axis
◆ SetFromAxisY()
|
inline |
Fill this 3x3 matrix as a rotation matrix with the Y axis along the provided direction.
Uses the Gram-Schmidt orthonormalization. The optional argument z_sugg, together with y_dir, suggests the YZ plane (as long as y_dir is not too close to z_sugg, in which case a different direction is selected).
- Parameters
-
y_dir Y axis z_sugg suggested Z axis
◆ SetFromAxisZ()
|
inline |
Fill this 3x3 matrix as a rotation matrix with the Z axis along the provided direction.
Uses the Gram-Schmidt orthonormalization. The optional argument x_sugg, together with z_dir, suggests the ZX plane (as long as z_dir is not too close to x_sugg, in which case a different direction is selected).
- Parameters
-
z_dir Z axis x_sugg suggested X axis
◆ SliceVector()
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.
◆ StreamOut()
|
inline |
Serialization of a sparse matrix to an ASCII stream (e.g., a file) in COO sparse matrix format.
By default, uses 0-based indices. If one_indexed=true, row and column indices start at 1 (as in Matlab).