![]() |
Eigen
3.4.90 (git rev 9589cc4e7fd8e4538bedef80dd36c7738977a8be)
|
This module provides support for:
Topics | |
Global aligned box typedefs | |
Classes | |
class | Eigen::AlignedBox< Scalar_, AmbientDim_ > |
An axis aligned box. More... | |
class | Eigen::AngleAxis< Scalar_ > |
Represents a 3D rotation as a rotation angle around an arbitrary 3D axis. More... | |
class | Eigen::Homogeneous< MatrixType, Direction_ > |
Expression of one (or a set of) homogeneous vector(s) More... | |
class | Eigen::Hyperplane< Scalar_, AmbientDim_, Options_ > |
A hyperplane. More... | |
class | Eigen::Map< const Quaternion< Scalar_ >, Options_ > |
Quaternion expression mapping a constant memory buffer. More... | |
class | Eigen::Map< Quaternion< Scalar_ >, Options_ > |
Expression of a quaternion from a memory buffer. More... | |
class | Eigen::ParametrizedLine< Scalar_, AmbientDim_, Options_ > |
A parametrized line. More... | |
class | Eigen::Quaternion< Scalar_, Options_ > |
The quaternion class used to represent 3D orientations and rotations. More... | |
class | Eigen::QuaternionBase< Derived > |
Base class for quaternion expressions. More... | |
class | Eigen::Rotation2D< Scalar_ > |
Represents a rotation/orientation in a 2 dimensional space. More... | |
class | Eigen::Transform< Scalar_, Dim_, Mode_, Options_ > |
Represents an homogeneous transformation in a N dimensional space. More... | |
class | Eigen::Translation< Scalar_, Dim_ > |
Represents a translation transformation. More... | |
class | Eigen::UniformScaling< Scalar_ > |
Represents a generic uniform scaling transformation. More... | |
Typedefs | |
typedef AngleAxis< double > | Eigen::AngleAxisd |
typedef AngleAxis< float > | Eigen::AngleAxisf |
typedef Quaternion< double > | Eigen::Quaterniond |
typedef Quaternion< float > | Eigen::Quaternionf |
typedef Map< Quaternion< double >, Aligned > | Eigen::QuaternionMapAlignedd |
typedef Map< Quaternion< float >, Aligned > | Eigen::QuaternionMapAlignedf |
typedef Map< Quaternion< double >, 0 > | Eigen::QuaternionMapd |
typedef Map< Quaternion< float >, 0 > | Eigen::QuaternionMapf |
typedef Rotation2D< double > | Eigen::Rotation2Dd |
typedef Rotation2D< float > | Eigen::Rotation2Df |
Functions | |
Matrix< Scalar, 3, 1 > | Eigen::MatrixBase< Derived >::canonicalEulerAngles (Index a0, Index a1, Index a2) const |
template<typename OtherDerived> | |
internal::cross_impl< Derived, OtherDerived >::return_type | Eigen::MatrixBase< Derived >::cross (const MatrixBase< OtherDerived > &other) const |
template<typename OtherDerived> | |
const CrossReturnType | Eigen::VectorwiseOp< ExpressionType, Direction >::cross (const MatrixBase< OtherDerived > &other) const |
template<typename OtherDerived> | |
PlainObject | Eigen::MatrixBase< Derived >::cross3 (const MatrixBase< OtherDerived > &other) const |
EIGEN_DEPRECATED Matrix< Scalar, 3, 1 > | Eigen::MatrixBase< Derived >::eulerAngles (Index a0, Index a1, Index a2) const |
const HNormalizedReturnType | Eigen::MatrixBase< Derived >::hnormalized () const |
homogeneous normalization | |
const HNormalizedReturnType | Eigen::VectorwiseOp< ExpressionType, Direction >::hnormalized () const |
column or row-wise homogeneous normalization | |
HomogeneousReturnType | Eigen::MatrixBase< Derived >::homogeneous () const |
HomogeneousReturnType | Eigen::VectorwiseOp< ExpressionType, Direction >::homogeneous () const |
template<typename Derived, typename OtherDerived> | |
internal::umeyama_transform_matrix_type< Derived, OtherDerived >::type | Eigen::umeyama (const MatrixBase< Derived > &src, const MatrixBase< OtherDerived > &dst, bool with_scaling=true) |
Returns the transformation between two point sets. | |
PlainObject | Eigen::MatrixBase< Derived >::unitOrthogonal (void) const |
typedef AngleAxis<double> Eigen::AngleAxisd |
double precision angle-axis type
typedef AngleAxis<float> Eigen::AngleAxisf |
single precision angle-axis type
typedef Quaternion<double> Eigen::Quaterniond |
double precision quaternion type
typedef Quaternion<float> Eigen::Quaternionf |
single precision quaternion type
typedef Map<Quaternion<double>, Aligned> Eigen::QuaternionMapAlignedd |
Map a 16-byte aligned array of double precision scalars as a quaternion
typedef Map<Quaternion<float>, Aligned> Eigen::QuaternionMapAlignedf |
Map a 16-byte aligned array of single precision scalars as a quaternion
typedef Map<Quaternion<double>, 0> Eigen::QuaternionMapd |
Map an unaligned array of double precision scalars as a quaternion
typedef Map<Quaternion<float>, 0> Eigen::QuaternionMapf |
Map an unaligned array of single precision scalars as a quaternion
typedef Rotation2D<double> Eigen::Rotation2Dd |
double precision 2D rotation type
typedef Rotation2D<float> Eigen::Rotation2Df |
single precision 2D rotation type
|
inline |
This is defined in the Geometry module.
*this
using the convention defined by the triplet (a0,a1,a2)Each of the three parameters a0,a1,a2 represents the respective rotation axis as an integer in {0,1,2}. For instance, in:
"2" represents the z axis and "0" the x axis, etc. The returned angles are such that we have the following equality:
This corresponds to the right-multiply conventions (with right hand side frames).
For Tait-Bryan angle configurations (a0 != a2), the returned angles are in the ranges [-pi:pi]x[-pi/2:pi/2]x[-pi:pi]. For proper Euler angle configurations (a0 == a2), the returned angles are in the ranges [-pi:pi]x[0:pi]x[-pi:pi].
The approach used is also described here: https://d3cw3dd2w32x2b.cloudfront.net/wp-content/uploads/2012/07/euler-angles.pdf
|
inline |
This is defined in the Geometry module.
*this
and other. This is either a scalar for size-2 vectors or a size-3 vector for size-3 vectors.This method is implemented for two different cases: between vectors of fixed size 2 and between vectors of fixed size 3.
For vectors of size 3, the output is simply the traditional cross product.
For vectors of size 2, the output is a scalar. Given vectors
const VectorwiseOp< ExpressionType, Direction >::CrossReturnType Eigen::VectorwiseOp< ExpressionType, Direction >::cross | ( | const MatrixBase< OtherDerived > & | other | ) | const |
This is defined in the Geometry module.
The referenced matrix must have one dimension equal to 3. The result matrix has the same dimensions than the referenced one.
|
inline |
This is defined in the Geometry module.
*this
and other using only the x, y, and z coefficientsThe size of *this
and other must be four. This function is especially useful when using 4D vectors instead of 3D ones to get advantage of SSE/AltiVec vectorization.
|
inline |
This is defined in the Geometry module.
*this
using the convention defined by the triplet (a0,a1,a2)NB: The returned angles are in non-canonical ranges [0:pi]x[-pi:pi]x[-pi:pi]. For canonical Tait-Bryan/proper Euler ranges, use canonicalEulerAngles.
|
inline |
homogeneous normalization
This is defined in the Geometry module.
*this
divided by that last coefficient.This can be used to convert homogeneous coordinates to affine coordinates.
It is essentially a shortcut for:
Example:
Output:
v = 0.696 0.205 -0.415 0.334]^T v.hnormalized() = 2.08 0.614 -1.24]^T P*v = -0.621 0.974 0.15 0.00432]^T (P*v).hnormalized() = -144 225 34.6]^T
|
inline |
column or row-wise homogeneous normalization
This is defined in the Geometry module.
*this
divided by the last coefficient of each column (or row).This can be used to convert homogeneous coordinates to affine coordinates.
It is conceptually equivalent to calling MatrixBase::hnormalized() to each column (or row) of *this
.
Example:
Output:
The matrix M is: 0.696 -0.47 0.0241 0.134 -0.727 0.205 0.928 0.0723 -0.16 0.74 -0.415 0.445 0.432 -0.00986 -0.757 0.334 -0.633 -0.046 -0.498 0.741 M.colwise().hnormalized(): 2.08 0.742 -0.524 -0.269 -0.982 0.614 -1.47 -1.57 0.32 0.999 -1.24 -0.704 -9.39 0.0198 -1.02 P*M: -0.17 -0.888 -0.0781 -0.119 -0.106 0.343 -0.998 0.0615 0.238 -1.36 0.246 -0.985 0.317 -0.0312 -1.41 -0.73 0.168 -0.165 0.242 0.0897 (P*M).colwise().hnormalized(): 0.233 -5.3 0.473 -0.491 -1.18 -0.47 -5.96 -0.372 0.986 -15.1 -0.336 -5.88 -1.92 -0.129 -15.8
|
inline |
This is defined in the Geometry module.
This can be used to convert affine coordinates to homogeneous coordinates.
This is only for vectors (either row-vectors or column-vectors), i.e. matrices which are known at compile-time to have either one row or one column.
Example:
Output:
v = [ 0.696 0.205 -0.415]^T h.homogeneous() = [ 0.696 0.205 -0.415 1]^T (P * v.homogeneous()) = [-0.376 -1.1 1.47 -0.354]^T (P * v.homogeneous()).hnormalized() = [ 1.06 3.12 -4.14]^T
|
inline |
This is defined in the Geometry module.
This can be used to convert affine coordinates to homogeneous coordinates.
Example:
Output:
The matrix M is: 0.696 0.334 0.445 0.0723 0.134 0.205 -0.47 -0.633 0.432 -0.16 -0.415 0.928 0.0241 -0.046 -0.00986 M.colwise().homogeneous(): 0.696 0.334 0.445 0.0723 0.134 0.205 -0.47 -0.633 0.432 -0.16 -0.415 0.928 0.0241 -0.046 -0.00986 1 1 1 1 1 P * M.colwise().homogeneous(): -0.316 -1.75 -1.18 -0.145 -0.644 -0.221 -0.856 -0.198 -0.103 -0.0481 1.22 -0.635 -0.00766 0.677 0.191 0.798 -0.446 -0.0232 1.2 0.631 P * M.colwise().homogeneous().hnormalized(): -0.396 3.94 50.8 -0.121 -1.02 -0.277 1.92 8.55 -0.0861 -0.0763 1.53 1.42 0.33 0.565 0.303
internal::umeyama_transform_matrix_type< Derived, OtherDerived >::type Eigen::umeyama | ( | const MatrixBase< Derived > & | src, |
const MatrixBase< OtherDerived > & | dst, | ||
bool | with_scaling = true ) |
Returns the transformation between two point sets.
This is defined in the Geometry module.
The algorithm is based on: "Least-squares estimation of transformation parameters between two point patterns", Shinji Umeyama, PAMI 1991, DOI: 10.1109/34.88573
It estimates parameters
is minimized.
The algorithm is based on the analysis of the covariance matrix
Currently the method is working only for floating point matrices.
src | Source points ![]() |
dst | Destination points ![]() |
with_scaling | Sets ![]() false is passed. |
|
inline |
This is defined in the Geometry module.
*this
The size of *this
must be at least 2. If the size is exactly 2, then the returned vector is a counter clock wise rotation of *this
, i.e., (-y,x).normalized().