Medial Code Documentation
Loading...
Searching...
No Matches
Modules | Data Structures | Typedefs | Functions
Geometry_Module

Modules

 Global aligned box typedefs
 Eigen defines several typedef shortcuts for most common aligned box types.
 

Data Structures

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::AlignedBox< _Scalar, _AmbientDim >
 \geometry_module More...
 
class  Eigen::AngleAxis< _Scalar >
 \geometry_module More...
 
class  Eigen::Homogeneous< MatrixType, _Direction >
 \geometry_module More...
 
class  Eigen::Hyperplane< _Scalar, _AmbientDim, _Options >
 \geometry_module More...
 
class  Eigen::ParametrizedLine< _Scalar, _AmbientDim, _Options >
 \geometry_module More...
 
class  Eigen::QuaternionBase< Derived >
 \geometry_module More...
 
class  Eigen::Quaternion< _Scalar, _Options >
 \geometry_module More...
 
class  Eigen::Rotation2D< _Scalar >
 \geometry_module More...
 
class  Eigen::UniformScaling< _Scalar >
 \geometry_module More...
 
class  Eigen::Transform< _Scalar, _Dim, _Mode, _Options >
 \geometry_module More...
 
class  Eigen::Translation< _Scalar, _Dim >
 \geometry_module More...
 

Typedefs

typedef AngleAxis< float > Eigen::AngleAxisf
 single precision angle-axis type
 
typedef AngleAxis< doubleEigen::AngleAxisd
 double precision angle-axis type
 
typedef Quaternion< float > Eigen::Quaternionf
 single precision quaternion type
 
typedef Quaternion< doubleEigen::Quaterniond
 double precision quaternion type
 
typedef Map< Quaternion< float >, 0 > Eigen::QuaternionMapf
 Map an unaligned 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 >, Aligned > Eigen::QuaternionMapAlignedf
 Map a 16-byte aligned array of single precision scalars as a quaternion.
 
typedef Map< Quaternion< double >, Aligned > Eigen::QuaternionMapAlignedd
 Map a 16-byte aligned array of double precision scalars as a quaternion.
 
typedef Rotation2D< float > Eigen::Rotation2Df
 single precision 2D rotation type
 
typedef Rotation2D< doubleEigen::Rotation2Dd
 double precision 2D rotation type
 
typedef Transform< float, 2, Isometry > Eigen::Isometry2f
 
typedef Transform< float, 3, Isometry > Eigen::Isometry3f
 
typedef Transform< double, 2, Isometry > Eigen::Isometry2d
 
typedef Transform< double, 3, Isometry > Eigen::Isometry3d
 
typedef Transform< float, 2, Affine > Eigen::Affine2f
 
typedef Transform< float, 3, Affine > Eigen::Affine3f
 
typedef Transform< double, 2, Affine > Eigen::Affine2d
 
typedef Transform< double, 3, Affine > Eigen::Affine3d
 
typedef Transform< float, 2, AffineCompact > Eigen::AffineCompact2f
 
typedef Transform< float, 3, AffineCompact > Eigen::AffineCompact3f
 
typedef Transform< double, 2, AffineCompact > Eigen::AffineCompact2d
 
typedef Transform< double, 3, AffineCompact > Eigen::AffineCompact3d
 
typedef Transform< float, 2, Projective > Eigen::Projective2f
 
typedef Transform< float, 3, Projective > Eigen::Projective3f
 
typedef Transform< double, 2, Projective > Eigen::Projective2d
 
typedef Transform< double, 3, Projective > Eigen::Projective3d
 

Functions

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)
 \geometry_module
 
EIGEN_DEVICE_FUNC Matrix< Scalar, 3, 1 > Eigen::MatrixBase< Derived >::eulerAngles (Index a0, Index a1, Index a2) const
 \geometry_module
 
EIGEN_DEVICE_FUNC HomogeneousReturnType Eigen::MatrixBase< Derived >::homogeneous () const
 \geometry_module
 
EIGEN_DEVICE_FUNC HomogeneousReturnType Eigen::VectorwiseOp< ExpressionType, Direction >::homogeneous () const
 \geometry_module
 
EIGEN_DEVICE_FUNC const HNormalizedReturnType Eigen::MatrixBase< Derived >::hnormalized () const
 \geometry_module
 
EIGEN_DEVICE_FUNC const HNormalizedReturnType Eigen::VectorwiseOp< ExpressionType, Direction >::hnormalized () const
 \geometry_module
 
template<typename OtherDerived >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE MatrixBase< Derived >::template cross_product_return_type< OtherDerived >::type Eigen::MatrixBase< Derived >::cross (const MatrixBase< OtherDerived > &other) const
 \geometry_module
 
template<typename OtherDerived >
EIGEN_DEVICE_FUNC PlainObject Eigen::MatrixBase< Derived >::cross3 (const MatrixBase< OtherDerived > &other) const
 \geometry_module
 
template<typename OtherDerived >
EIGEN_DEVICE_FUNC const CrossReturnType Eigen::VectorwiseOp< ExpressionType, Direction >::cross (const MatrixBase< OtherDerived > &other) const
 \geometry_module
 
EIGEN_DEVICE_FUNC PlainObject Eigen::MatrixBase< Derived >::unitOrthogonal (void) const
 \geometry_module
 

Detailed Description

Function Documentation

◆ cross() [1/2]

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

\geometry_module

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()

◆ cross() [2/2]

template<typename Derived >
template<typename OtherDerived >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE MatrixBase< Derived >::template cross_product_return_type< OtherDerived >::type Eigen::MatrixBase< Derived >::cross ( const MatrixBase< OtherDerived > &  other) const

\geometry_module

Returns
the cross product of *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 $ (\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()

◆ cross3()

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

\geometry_module

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_DEVICE_FUNC Matrix< typename MatrixBase< Derived >::Scalar, 3, 1 > Eigen::MatrixBase< Derived >::eulerAngles ( Index  a0,
Index  a1,
Index  a2 
) const
inline

\geometry_module

Returns
the 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);
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
EIGEN_DEVICE_FUNC Matrix< Scalar, 3, 1 > eulerAngles(Index a0, Index a1, Index a2) const
\geometry_module
Definition EulerAngles.h:37

"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())
* AngleAxisf(ea[1], Vector3f::UnitX())
* AngleAxisf(ea[2], Vector3f::UnitZ());
AngleAxis< float > AngleAxisf
single precision angle-axis type
Definition AngleAxis.h:157

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

See also
class AngleAxis

◆ hnormalized() [1/2]

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

\geometry_module

homogeneous normalization

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);

Example:

Output:

See also
VectorwiseOp::hnormalized()

◆ hnormalized() [2/2]

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

\geometry_module

column or row-wise homogeneous normalization

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:

Output:

See also
MatrixBase::hnormalized()

◆ homogeneous() [1/2]

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

\geometry_module

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.

\only_for_vectors

Example:

Output:

See also
VectorwiseOp::homogeneous(), class Homogeneous

◆ homogeneous() [2/2]

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

\geometry_module

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:

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 
)

\geometry_module

Returns the transformation between two point sets.

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.

Todo:
Should the return type of umeyama() become a Transform?
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 >
EIGEN_DEVICE_FUNC MatrixBase< Derived >::PlainObject Eigen::MatrixBase< Derived >::unitOrthogonal ( void  ) const
inline

\geometry_module

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()