Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
Geometry module

Detailed Description

This module provides support for:

#include <Eigen/Geometry>

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 >
std::conditional_t< SizeAtCompileTime==2, Scalar, 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
 
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 Documentation

◆ AngleAxisd

typedef AngleAxis<double> Eigen::AngleAxisd

double precision angle-axis type

◆ AngleAxisf

typedef AngleAxis<float> Eigen::AngleAxisf

single precision angle-axis type

◆ Quaterniond

typedef Quaternion<double> Eigen::Quaterniond

double precision quaternion type

◆ Quaternionf

single precision quaternion type

◆ QuaternionMapAlignedd

typedef Map<Quaternion<double>, Aligned> Eigen::QuaternionMapAlignedd

Map a 16-byte aligned array of double precision scalars as a quaternion

◆ QuaternionMapAlignedf

typedef Map<Quaternion<float>, Aligned> Eigen::QuaternionMapAlignedf

Map a 16-byte aligned array of single precision scalars as a quaternion

◆ QuaternionMapd

typedef Map<Quaternion<double>, 0> Eigen::QuaternionMapd

Map an unaligned array of double precision scalars as a quaternion

◆ QuaternionMapf

typedef Map<Quaternion<float>, 0> Eigen::QuaternionMapf

Map an unaligned array of single precision scalars as a quaternion

◆ Rotation2Dd

typedef Rotation2D<double> Eigen::Rotation2Dd

double precision 2D rotation type

◆ Rotation2Df

single precision 2D rotation type

Function Documentation

◆ canonicalEulerAngles()

template<typename Derived >
Matrix< typename MatrixBase< Derived >::Scalar, 3, 1 > Eigen::MatrixBase< Derived >::canonicalEulerAngles ( Index a0,
Index a1,
Index a2 ) const
inline

This is defined in the Geometry module.

#include <Eigen/Geometry>
Returns
the canonical Euler-angles of the rotation matrix *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:

Vector3f ea = mat.eulerAngles(2, 0, 2);
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
EIGEN_DEPRECATED Matrix< Scalar, 3, 1 > eulerAngles(Index a0, Index a1, Index a2) const
Definition EulerAngles.h:137

"2" represents the z axis and "0" the x axis, etc. The returned angles are such that we have the following equality:

mat == AngleAxisf(ea[0], Vector3f::UnitZ())
static const BasisReturnType UnitX()
Definition CwiseNullaryOp.h:895
static const BasisReturnType UnitZ()
Definition CwiseNullaryOp.h:919
AngleAxis< float > AngleAxisf
Definition AngleAxis.h:165

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

See also
class AngleAxis

◆ cross() [1/2]

template<typename Derived >
template<typename OtherDerived >
std::conditional_t< SizeAtCompileTime==2, Scalar, PlainObject > Eigen::MatrixBase< Derived >::cross ( const MatrixBase< OtherDerived > & other) const
inline

This is defined in the Geometry module.

#include <Eigen/Geometry>
Returns
the cross product of *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 \( v = \begin{bmatrix} v_1 & v_2 \end{bmatrix} \) and \( w = \begin{bmatrix} w_1 & w_2 \end{bmatrix} \), the result is simply \( v\times w = \overline{v_1 w_2 - v_2 w_1} = \text{conj}\left|\begin{smallmatrix} v_1 & w_1 \\ v_2 & w_2 \end{smallmatrix}\right| \); or, to put it differently, it is the third coordinate of the cross product of \( \begin{bmatrix} v_1 & v_2 & v_3 \end{bmatrix} \) and \( \begin{bmatrix} w_1 & w_2 & w_3 \end{bmatrix} \). For real-valued inputs, the result can be interpreted as the signed area of a parallelogram spanned by the two vectors.

Note
With complex numbers, the cross product is implemented as \( (\mathbf{a}+i\mathbf{b}) \times (\mathbf{c}+i\mathbf{d}) = (\mathbf{a} \times \mathbf{c} - \mathbf{b} \times \mathbf{d}) - i(\mathbf{a} \times \mathbf{d} + \mathbf{b} \times \mathbf{c})\)
See also
MatrixBase::cross3()

◆ cross() [2/2]

template<typename ExpressionType , int Direction>
template<typename OtherDerived >
const VectorwiseOp< ExpressionType, Direction >::CrossReturnType Eigen::VectorwiseOp< ExpressionType, Direction >::cross ( const MatrixBase< OtherDerived > & other) const

This is defined in the Geometry module.

#include <Eigen/Geometry>
Returns
a matrix expression of the cross product of each column or row of the referenced expression with the other vector.

The referenced matrix must have one dimension equal to 3. The result matrix has the same dimensions than the referenced one.

See also
MatrixBase::cross()

◆ cross3()

template<typename Derived >
template<typename OtherDerived >
MatrixBase< Derived >::PlainObject Eigen::MatrixBase< Derived >::cross3 ( const MatrixBase< OtherDerived > & other) const
inline

This is defined in the Geometry module.

#include <Eigen/Geometry>
Returns
the cross product of *this and other using only the x, y, and z coefficients

The 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.

See also
MatrixBase::cross()

◆ eulerAngles()

template<typename Derived >
EIGEN_DEPRECATED Matrix< typename MatrixBase< Derived >::Scalar, 3, 1 > Eigen::MatrixBase< Derived >::eulerAngles ( Index a0,
Index a1,
Index a2 ) const
inline

This is defined in the Geometry module.

#include <Eigen/Geometry>
Returns
the Euler-angles of the rotation matrix *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.

See also
MatrixBase::canonicalEulerAngles
class AngleAxis

◆ hnormalized() [1/2]

template<typename Derived >
const MatrixBase< Derived >::HNormalizedReturnType Eigen::MatrixBase< Derived >::hnormalized ( ) const
inline

homogeneous normalization

This is defined in the Geometry module.

#include <Eigen/Geometry>
Returns
a vector expression of the N-1 first coefficients of *this divided by that last coefficient.

This can be used to convert homogeneous coordinates to affine coordinates.

It is essentially a shortcut for:

this->head(this->size()-1)/this->coeff(this->size()-1);
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition EigenBase.h:64
EIGEN_CONSTEXPR CoeffReturnType coeff(Index row, Index col) const
Definition DenseCoeffsBase.h:92

Example:

Vector4d v = Vector4d::Random();
Projective3d P(Matrix4d::Random());
cout << "v = " << v.transpose() << "]^T" << endl;
cout << "v.hnormalized() = " << v.hnormalized().transpose() << "]^T" << endl;
cout << "P*v = " << (P * v).transpose() << "]^T" << endl;
cout << "(P*v).hnormalized() = " << (P * v).hnormalized().transpose() << "]^T" << endl;
TransposeReturnType transpose()
Definition Transpose.h:160
Represents an homogeneous transformation in a N dimensional space.
Definition Transform.h:192
const HNormalizedReturnType hnormalized() const
homogeneous normalization
Definition Homogeneous.h:166

Output:

See also
VectorwiseOp::hnormalized()

◆ hnormalized() [2/2]

template<typename ExpressionType , int Direction>
const VectorwiseOp< ExpressionType, Direction >::HNormalizedReturnType Eigen::VectorwiseOp< ExpressionType, Direction >::hnormalized ( ) const
inline

column or row-wise homogeneous normalization

This is defined in the Geometry module.

#include <Eigen/Geometry>
Returns
an expression of the first N-1 coefficients of each column (or row) of *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:

Matrix4Xd M = Matrix4Xd::Random(4, 5);
Projective3d P(Matrix4d::Random());
cout << "The matrix M is:" << endl << M << endl << endl;
cout << "M.colwise().hnormalized():" << endl << M.colwise().hnormalized() << endl << endl;
cout << "P*M:" << endl << P * M << endl << endl;
cout << "(P*M).colwise().hnormalized():" << endl << (P * M).colwise().hnormalized() << endl << endl;
ConstColwiseReturnType colwise() const
Definition DenseBase.h:511
const HNormalizedReturnType hnormalized() const
column or row-wise homogeneous normalization
Definition Homogeneous.h:191

Output:

See also
MatrixBase::hnormalized()

◆ homogeneous() [1/2]

template<typename Derived >
MatrixBase< Derived >::HomogeneousReturnType Eigen::MatrixBase< Derived >::homogeneous ( ) const
inline

This is defined in the Geometry module.

#include <Eigen/Geometry>
Returns
a vector expression that is one longer than the vector argument, with the value 1 symbolically appended as the last coefficient.

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:

Vector3d v = Vector3d::Random(), w;
Projective3d P(Matrix4d::Random());
cout << "v = [" << v.transpose() << "]^T" << endl;
cout << "h.homogeneous() = [" << v.homogeneous().transpose() << "]^T" << endl;
cout << "(P * v.homogeneous()) = [" << (P * v.homogeneous()).transpose() << "]^T" << endl;
cout << "(P * v.homogeneous()).hnormalized() = [" << (P * v.homogeneous()).eval().hnormalized().transpose() << "]^T"
<< endl;
HomogeneousReturnType homogeneous() const
Definition Homogeneous.h:126

Output:

See also
VectorwiseOp::homogeneous(), class Homogeneous

◆ homogeneous() [2/2]

template<typename ExpressionType , int Direction>
Homogeneous< ExpressionType, Direction > Eigen::VectorwiseOp< ExpressionType, Direction >::homogeneous ( ) const
inline

This is defined in the Geometry module.

#include <Eigen/Geometry>
Returns
an expression where the value 1 is symbolically appended as the final coefficient to each column (or row) of the matrix.

This can be used to convert affine coordinates to homogeneous coordinates.

Example:

Matrix3Xd M = Matrix3Xd::Random(3, 5);
Projective3d P(Matrix4d::Random());
cout << "The matrix M is:" << endl << M << endl << endl;
cout << "M.colwise().homogeneous():" << endl << M.colwise().homogeneous() << endl << endl;
cout << "P * M.colwise().homogeneous():" << endl << P * M.colwise().homogeneous() << endl << endl;
cout << "P * M.colwise().homogeneous().hnormalized(): " << endl
<< (P * M.colwise().homogeneous()).colwise().hnormalized() << endl
<< endl;
HomogeneousReturnType homogeneous() const
Definition Homogeneous.h:143

Output:

See also
MatrixBase::homogeneous(), class Homogeneous

◆ umeyama()

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.

This is defined in the Geometry module.

#include <Eigen/Geometry>

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 \( c, \mathbf{R}, \) and \( \mathbf{t} \) such that

\begin{align*} \frac{1}{n} \sum_{i=1}^n \vert\vert y_i - (c\mathbf{R}x_i + \mathbf{t}) \vert\vert_2^2 \end{align*}

is minimized.

The algorithm is based on the analysis of the covariance matrix \( \Sigma_{\mathbf{x}\mathbf{y}} \in \mathbb{R}^{d \times d} \) of the input point sets \( \mathbf{x} \) and \( \mathbf{y} \) where \(d\) is corresponding to the dimension (which is typically small). The analysis is involving the SVD having a complexity of \(O(d^3)\) though the actual computational effort lies in the covariance matrix computation which has an asymptotic lower bound of \(O(dm)\) when the input point sets have dimension \(d \times m\).

Currently the method is working only for floating point matrices.

Parameters
srcSource points \( \mathbf{x} = \left( x_1, \hdots, x_n \right) \).
dstDestination points \( \mathbf{y} = \left( y_1, \hdots, y_n \right) \).
with_scalingSets \( c=1 \) when false is passed.
Returns
The homogeneous transformation

\begin{align*} T = \begin{bmatrix} c\mathbf{R} & \mathbf{t} \\ \mathbf{0} & 1 \end{bmatrix} \end{align*}

minimizing the residual above. This transformation is always returned as an Eigen::Matrix.

◆ unitOrthogonal()

template<typename Derived >
MatrixBase< Derived >::PlainObject Eigen::MatrixBase< Derived >::unitOrthogonal ( void ) const
inline

This is defined in the Geometry module.

#include <Eigen/Geometry>
Returns
a unit vector which is orthogonal to *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().

See also
cross()