Modules | Classes | Typedefs | Functions
Geometry_Module
Collaboration diagram for Geometry_Module:

Modules

 Global aligned box typedefs
 

Classes

class  Eigen::AlignedBox
 An axis aligned box. More...
 
class  Eigen::AngleAxis
 Represents a 3D rotation as a rotation angle around an arbitrary 3D axis. More...
 
class  Eigen::Homogeneous
 Expression of one (or a set of) homogeneous vector(s) More...
 
class  Eigen::Hyperplane
 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
 A parametrized line. More...
 
class  Eigen::Quaternion
 The quaternion class used to represent 3D orientations and rotations. More...
 
class  Eigen::QuaternionBase
 Base class for quaternion expressions. More...
 
class  Eigen::Rotation2D
 Represents a rotation/orientation in a 2 dimensional space. More...
 
class  Scaling
 Represents a generic uniform scaling transformation. More...
 
class  Eigen::Transform
 Represents an homogeneous transformation in a N dimensional space. More...
 
class  Eigen::Translation
 Represents a translation transformation. More...
 

Typedefs

typedef Transform< double, 2, AffineEigen::Affine2d
 
typedef Transform< float, 2, AffineEigen::Affine2f
 
typedef Transform< double, 3, AffineEigen::Affine3d
 
typedef Transform< float, 3, AffineEigen::Affine3f
 
typedef Transform< double, 2, AffineCompactEigen::AffineCompact2d
 
typedef Transform< float, 2, AffineCompactEigen::AffineCompact2f
 
typedef Transform< double, 3, AffineCompactEigen::AffineCompact3d
 
typedef Transform< float, 3, AffineCompactEigen::AffineCompact3f
 
typedef AngleAxis< double > Eigen::AngleAxisd
 
typedef AngleAxis< float > Eigen::AngleAxisf
 
typedef Transform< double, 2, IsometryEigen::Isometry2d
 
typedef Transform< float, 2, IsometryEigen::Isometry2f
 
typedef Transform< double, 3, IsometryEigen::Isometry3d
 
typedef Transform< float, 3, IsometryEigen::Isometry3f
 
typedef Transform< double, 2, ProjectiveEigen::Projective2d
 
typedef Transform< float, 2, ProjectiveEigen::Projective2f
 
typedef Transform< double, 3, ProjectiveEigen::Projective3d
 
typedef Transform< float, 3, ProjectiveEigen::Projective3f
 
typedef Quaternion< double > Eigen::Quaterniond
 
typedef Quaternion< float > Eigen::Quaternionf
 
typedef Map< Quaternion< double >, AlignedEigen::QuaternionMapAlignedd
 
typedef Map< Quaternion< float >, AlignedEigen::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 >
EIGEN_DEVICE_FUNC MatrixBase< Derived >::template cross_product_return_type< OtherDerived >::type Eigen::MatrixBase< Derived >::cross (const MatrixBase< OtherDerived > &other) const
 
template<typename OtherDerived >
const EIGEN_DEVICE_FUNC CrossReturnType Eigen::VectorwiseOp::cross (const MatrixBase< OtherDerived > &other) const
 
template<typename OtherDerived >
EIGEN_DEVICE_FUNC PlainObject Eigen::MatrixBase< Derived >::cross3 (const MatrixBase< OtherDerived > &other) const
 
EIGEN_DEVICE_FUNC Matrix< Scalar, 3, 1 > Eigen::MatrixBase< Derived >::eulerAngles (Index a0, Index a1, Index a2) const
 
const EIGEN_DEVICE_FUNC HNormalizedReturnType Eigen::MatrixBase< Derived >::hnormalized () const
 homogeneous normalization More...
 
const EIGEN_DEVICE_FUNC HNormalizedReturnType Eigen::VectorwiseOp::hnormalized () const
 column or row-wise homogeneous normalization More...
 
EIGEN_DEVICE_FUNC HomogeneousReturnType Eigen::MatrixBase< Derived >::homogeneous () const
 
EIGEN_DEVICE_FUNC HomogeneousReturnType Eigen::VectorwiseOp::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. More...
 
EIGEN_DEVICE_FUNC PlainObject Eigen::MatrixBase< Derived >::unitOrthogonal (void) const
 

Detailed Description

Typedef Documentation

◆ Affine2d

typedef Transform<double,2,Affine> Eigen::Affine2d

Definition at line 712 of file Transform.h.

◆ Affine2f

typedef Transform<float,2,Affine> Eigen::Affine2f

Definition at line 708 of file Transform.h.

◆ Affine3d

typedef Transform<double,3,Affine> Eigen::Affine3d

Definition at line 714 of file Transform.h.

◆ Affine3f

typedef Transform<float,3,Affine> Eigen::Affine3f

Definition at line 710 of file Transform.h.

◆ AffineCompact2d

Definition at line 721 of file Transform.h.

◆ AffineCompact2f

Definition at line 717 of file Transform.h.

◆ AffineCompact3d

Definition at line 723 of file Transform.h.

◆ AffineCompact3f

Definition at line 719 of file Transform.h.

◆ AngleAxisd

typedef AngleAxis<double> Eigen::AngleAxisd

double precision angle-axis type

Definition at line 160 of file AngleAxis.h.

◆ AngleAxisf

typedef AngleAxis<float> Eigen::AngleAxisf

single precision angle-axis type

Definition at line 157 of file AngleAxis.h.

◆ Isometry2d

Definition at line 703 of file Transform.h.

◆ Isometry2f

Definition at line 699 of file Transform.h.

◆ Isometry3d

Definition at line 705 of file Transform.h.

◆ Isometry3f

Definition at line 701 of file Transform.h.

◆ Projective2d

Definition at line 730 of file Transform.h.

◆ Projective2f

Definition at line 726 of file Transform.h.

◆ Projective3d

Definition at line 732 of file Transform.h.

◆ Projective3f

Definition at line 728 of file Transform.h.

◆ Quaterniond

typedef Quaternion<double> Eigen::Quaterniond

double precision quaternion type

Definition at line 310 of file Quaternion.h.

◆ Quaternionf

single precision quaternion type

Definition at line 307 of file Quaternion.h.

◆ QuaternionMapAlignedd

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

Definition at line 422 of file Quaternion.h.

◆ QuaternionMapAlignedf

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

Definition at line 419 of file Quaternion.h.

◆ QuaternionMapd

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

Map an unaligned array of double precision scalars as a quaternion

Definition at line 416 of file Quaternion.h.

◆ QuaternionMapf

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

Map an unaligned array of single precision scalars as a quaternion

Definition at line 413 of file Quaternion.h.

◆ Rotation2Dd

typedef Rotation2D<double> Eigen::Rotation2Dd

double precision 2D rotation type

Definition at line 168 of file Rotation2D.h.

◆ Rotation2Df

single precision 2D rotation type

Definition at line 165 of file Rotation2D.h.

Function Documentation

◆ cross() [1/2]

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

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

Definition at line 34 of file OrthoMethods.h.

◆ cross() [2/2]

template<typename OtherDerived >
const EIGEN_DEVICE_FUNC VectorwiseOp< ExpressionType, Direction >::CrossReturnType Eigen::VectorwiseOp::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()

Definition at line 110 of file OrthoMethods.h.

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

Definition at line 82 of file OrthoMethods.h.

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

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

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

Definition at line 37 of file Eigen/src/Geometry/EulerAngles.h.

◆ hnormalized() [1/2]

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

homogeneous normalization

\geometry_module

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

Definition at line 172 of file Homogeneous.h.

◆ hnormalized() [2/2]

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

column or row-wise homogeneous normalization

\geometry_module

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

Definition at line 196 of file Homogeneous.h.

◆ homogeneous() [1/2]

template<typename Derived >
EIGEN_DEVICE_FUNC MatrixBase< Derived >::HomogeneousReturnType Eigen::MatrixBase< Derived >::homogeneous
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

Definition at line 130 of file Homogeneous.h.

◆ homogeneous() [2/2]

EIGEN_DEVICE_FUNC Homogeneous< ExpressionType, Direction > Eigen::VectorwiseOp::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

Definition at line 148 of file Homogeneous.h.

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

\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 $ 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 resudiual above. This transformation is always returned as an Eigen::Matrix.

Definition at line 95 of file Umeyama.h.

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

Definition at line 226 of file OrthoMethods.h.

head
EIGEN_DEVICE_FUNC SegmentReturnType head(Index n)
This is the const version of head(Index).
Definition: BlockMethods.h:919
Eigen::AngleAxisf
AngleAxis< float > AngleAxisf
Definition: AngleAxis.h:157
utility::tuple::size
static constexpr size_t size(Tuple< Args... > &)
Provides access to the number of elements in a tuple as a compile-time constant expression.
Definition: TensorSyclTuple.h:143
mat
else mat
Definition: eigenvalues.cpp:43


control_box_rst
Author(s): Christoph Rösmann
autogenerated on Wed Mar 2 2022 00:07:18