![]() |
Eigen
3.3.0
|
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 | Scaling |
Represents a generic uniform scaling transformation. 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... | |
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 | |
template<typename OtherDerived > | |
MatrixBase< Derived >::PlainObject | Eigen::MatrixBase< Derived >::cross (const MatrixBase< OtherDerived > &other) const |
template<typename OtherDerived > | |
const VectorwiseOp< ExpressionType, Direction >::CrossReturnType | Eigen::VectorwiseOp< ExpressionType, Direction >::cross (const MatrixBase< OtherDerived > &other) const |
template<typename OtherDerived > | |
MatrixBase< Derived >::PlainObject | Eigen::MatrixBase< Derived >::cross3 (const MatrixBase< OtherDerived > &other) const |
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
and other Here is a very good explanation of cross-product: http://xkcd.com/199/
With complex numbers, the cross product is implemented as
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)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).
The returned angles are in the ranges [0:pi]x[-pi:pi]x[-pi:pi].
|
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.0277 -0.649 -0.383 0.0691]^T v.hnormalized() = 0.402 -9.39 -5.54]^T P*v = 0.239 0.59 0.218 0.14]^T (P*v).hnormalized() = 1.71 4.23 1.56]^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.0277 0.895 -0.0105 -0.446 0.531 -0.649 -0.657 -0.751 -0.264 0.293 -0.383 0.404 -0.832 0.967 0.534 0.0691 -0.547 -0.221 0.0708 0.56 M.colwise().hnormalized(): 0.402 -1.64 0.0474 -6.29 0.948 -9.39 1.2 3.4 -3.73 0.523 -5.54 -0.739 3.77 13.7 0.953 P*M: 0.185 0.465 -0.116 0.0514 0.778 -0.755 -1.13 -1.17 0.667 0.375 -0.313 0.265 -0.865 0.749 0.896 -0.154 0.00761 -0.5 0.944 0.232 (P*M).colwise().hnormalized(): -1.2 61.1 0.233 0.0544 3.35 4.9 -149 2.34 0.706 1.62 2.04 34.8 1.73 0.793 3.86
|
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.0277 -0.649 -0.383]^T h.homogeneous() = [0.0277 -0.649 -0.383 1]^T (P * v.homogeneous()) = [0.512 0.733 0.862 0.715]^T (P * v.homogeneous()).hnormalized() = [0.716 1.03 1.21]^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.0277 0.0691 0.404 -0.751 -0.446 -0.649 0.895 -0.547 -0.832 -0.264 -0.383 -0.657 -0.0105 -0.221 0.967 M.colwise().homogeneous(): 0.0277 0.0691 0.404 -0.751 -0.446 -0.649 0.895 -0.547 -0.832 -0.264 -0.383 -0.657 -0.0105 -0.221 0.967 1 1 1 1 1 P * M.colwise().homogeneous(): -0.417 0.552 -0.472 -0.635 -0.735 -0.0734 1.03 0.0782 -0.655 -0.49 0.711 -0.58 1.06 0.746 1.43 0.7 1.1 0.941 0.244 0.597 P * M.colwise().homogeneous().hnormalized(): -0.596 0.503 -0.501 -2.6 -1.23 -0.105 0.937 0.0831 -2.68 -0.82 1.02 -0.528 1.13 3.05 2.4
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().