Geometry module

Modules

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

Matrix< Scalar, 3, 1 > Eigen::MatrixBase< Derived >::canonicalEulerAngles (Index a0, Index a1, Index a2) const
 
template<typename OtherDerived >
const CrossReturnType Eigen::VectorwiseOp< ExpressionType, Direction >::cross (const MatrixBase< OtherDerived > &other) const
 
template<typename OtherDerived >
std::conditional_t< SizeAtCompileTime==2, Scalar, PlainObjectEigen::MatrixBase< Derived >::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 More...
 
const HNormalizedReturnType Eigen::VectorwiseOp< ExpressionType, Direction >::hnormalized () const
 column or row-wise homogeneous normalization More...
 
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. More...
 
PlainObject Eigen::MatrixBase< Derived >::unitOrthogonal (void) const
 

Detailed Description

This module provides support for:

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 162 of file AngleAxis.h.

◆ AngleAxisf

typedef AngleAxis<float> Eigen::AngleAxisf

single precision angle-axis type

Definition at line 159 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 372 of file Quaternion.h.

◆ Quaternionf

single precision quaternion type

Definition at line 369 of file Quaternion.h.

◆ QuaternionMapAlignedd

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

Definition at line 484 of file Quaternion.h.

◆ QuaternionMapAlignedf

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

Definition at line 481 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 478 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 475 of file Quaternion.h.

◆ Rotation2Dd

typedef Rotation2D<double> Eigen::Rotation2Dd

double precision 2D rotation type

Definition at line 170 of file Rotation2D.h.

◆ Rotation2Df

single precision 2D rotation type

Definition at line 167 of file Rotation2D.h.

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);
Matrix< float, 3, 1 > Vector3f
3×1 vector of type float.
Definition: Matrix.h:501

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

static const BasisReturnType UnitX()
static const BasisReturnType UnitZ()
AngleAxis< float > AngleAxisf
Definition: AngleAxis.h:159

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

Definition at line 43 of file EulerAngles.h.

44 {
45  /* Implemented from Graphics Gems IV */
47 
48  Matrix<Scalar, 3, 1> res;
49 
50  const Index odd = ((a0 + 1) % 3 == a1) ? 0 : 1;
51  const Index i = a0;
52  const Index j = (a0 + 1 + odd) % 3;
53  const Index k = (a0 + 2 - odd) % 3;
54 
55  if (a0 == a2)
56  {
57  // Proper Euler angles (same first and last axis).
58  // The i, j, k indices enable addressing the input matrix as the XYX archetype matrix (see Graphics Gems IV),
59  // where e.g. coeff(k, i) means third column, first row in the XYX archetype matrix:
60  // c2 s2s1 s2c1
61  // s2s3 -c2s1s3 + c1c3 -c2c1s3 - s1c3
62  // -s2c3 c2s1c3 + c1s3 c2c1c3 - s1s3
63 
64  // Note: s2 is always positive.
65  Scalar s2 = numext::hypot(coeff(j, i), coeff(k, i));
66  if (odd)
67  {
68  res[0] = numext::atan2(coeff(j, i), coeff(k, i));
69  // s2 is always positive, so res[1] will be within the canonical [0, pi] range
70  res[1] = numext::atan2(s2, coeff(i, i));
71  }
72  else
73  {
74  // In the !odd case, signs of all three angles are flipped at the very end. To keep the solution within the canonical range,
75  // we flip the solution and make res[1] always negative here (since s2 is always positive, -atan2(s2, c2) will always be negative).
76  // The final flip at the end due to !odd will thus make res[1] positive and canonical.
77  // NB: in the general case, there are two correct solutions, but only one is canonical. For proper Euler angles,
78  // flipping from one solution to the other involves flipping the sign of the second angle res[1] and adding/subtracting pi
79  // to the first and third angles. The addition/subtraction of pi to the first angle res[0] is handled here by flipping
80  // the signs of arguments to atan2, while the calculation of the third angle does not need special adjustment since
81  // it uses the adjusted res[0] as the input and produces a correct result.
82  res[0] = numext::atan2(-coeff(j, i), -coeff(k, i));
83  res[1] = -numext::atan2(s2, coeff(i, i));
84  }
85 
86  // With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles,
87  // we can compute their respective rotation, and apply its inverse to M. Since the result must
88  // be a rotation around x, we have:
89  //
90  // c2 s1.s2 c1.s2 1 0 0
91  // 0 c1 -s1 * M = 0 c3 s3
92  // -s2 s1.c2 c1.c2 0 -s3 c3
93  //
94  // Thus: m11.c1 - m21.s1 = c3 & m12.c1 - m22.s1 = s3
95 
96  Scalar s1 = numext::sin(res[0]);
97  Scalar c1 = numext::cos(res[0]);
98  res[2] = numext::atan2(c1 * coeff(j, k) - s1 * coeff(k, k), c1 * coeff(j, j) - s1 * coeff(k, j));
99  }
100  else
101  {
102  // Tait-Bryan angles (all three axes are different; typically used for yaw-pitch-roll calculations).
103  // The i, j, k indices enable addressing the input matrix as the XYZ archetype matrix (see Graphics Gems IV),
104  // where e.g. coeff(k, i) means third column, first row in the XYZ archetype matrix:
105  // c2c3 s2s1c3 - c1s3 s2c1c3 + s1s3
106  // c2s3 s2s1s3 + c1c3 s2c1s3 - s1c3
107  // -s2 c2s1 c2c1
108 
109  res[0] = numext::atan2(coeff(j, k), coeff(k, k));
110 
111  Scalar c2 = numext::hypot(coeff(i, i), coeff(i, j));
112  // c2 is always positive, so the following atan2 will always return a result in the correct canonical middle angle range [-pi/2, pi/2]
113  res[1] = numext::atan2(-coeff(i, k), c2);
114 
115  Scalar s1 = numext::sin(res[0]);
116  Scalar c1 = numext::cos(res[0]);
117  res[2] = numext::atan2(s1 * coeff(k, i) - c1 * coeff(j, i), c1 * coeff(j, j) - s1 * coeff(k, j));
118  }
119  if (!odd)
120  {
121  res = -res;
122  }
123 
124  return res;
125 }
if((m *x).isApprox(y))
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
#define EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(TYPE, ROWS, COLS)
Definition: StaticAssert.h:56
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:61
CoeffReturnType coeff(Index row, Index col) const
EIGEN_ALWAYS_INLINE T sin(const T &x)
EIGEN_ALWAYS_INLINE T atan2(const T &y, const T &x)
EIGEN_ALWAYS_INLINE T cos(const T &x)
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
std::ptrdiff_t j

◆ cross() [1/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()

Definition at line 157 of file OrthoMethods.h.

158 {
160  EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename OtherDerived::Scalar>::value),
161  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
162 
163  typename internal::nested_eval<ExpressionType,2>::type mat(_expression());
164  typename internal::nested_eval<OtherDerived,2>::type vec(other.derived());
165 
167  if(Direction==Vertical)
168  {
169  eigen_assert(CrossReturnType::RowsAtCompileTime==3 && "the matrix must have exactly 3 rows");
170  res.row(0) = (mat.row(1) * vec.coeff(2) - mat.row(2) * vec.coeff(1)).conjugate();
171  res.row(1) = (mat.row(2) * vec.coeff(0) - mat.row(0) * vec.coeff(2)).conjugate();
172  res.row(2) = (mat.row(0) * vec.coeff(1) - mat.row(1) * vec.coeff(0)).conjugate();
173  }
174  else
175  {
176  eigen_assert(CrossReturnType::ColsAtCompileTime==3 && "the matrix must have exactly 3 columns");
177  res.col(0) = (mat.col(1) * vec.coeff(2) - mat.col(2) * vec.coeff(1)).conjugate();
178  res.col(1) = (mat.col(2) * vec.coeff(0) - mat.col(0) * vec.coeff(2)).conjugate();
179  res.col(2) = (mat.col(0) * vec.coeff(1) - mat.col(1) * vec.coeff(0)).conjugate();
180  }
181  return res;
182 }
ConjugateReturnType conjugate() const
#define eigen_assert(x)
Definition: Macros.h:902
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26
#define EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(TYPE, SIZE)
Definition: StaticAssert.h:51
ExpressionType::PlainObject CrossReturnType
Definition: VectorwiseOp.h:713
ExpressionType::Scalar Scalar
Definition: VectorwiseOp.h:191
const ExpressionType & _expression() const
Definition: VectorwiseOp.h:273
@ Vertical
Definition: Constants.h:266

◆ cross() [2/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()

Definition at line 92 of file OrthoMethods.h.

93 {
94  return internal::cross_impl<Derived, OtherDerived>::run(*this, other);
95 }

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

Definition at line 129 of file OrthoMethods.h.

130 {
133 
134  typedef typename internal::nested_eval<Derived,2>::type DerivedNested;
135  typedef typename internal::nested_eval<OtherDerived,2>::type OtherDerivedNested;
136  DerivedNested lhs(derived());
137  OtherDerivedNested rhs(other.derived());
138 
139  return internal::cross3_impl<Architecture::Target,
140  internal::remove_all_t<DerivedNested>,
141  internal::remove_all_t<OtherDerivedNested>>::run(lhs,rhs);
142 }
typename remove_all< T >::type remove_all_t
Definition: Meta.h:119

◆ 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

Definition at line 139 of file EulerAngles.h.

140 {
141  /* Implemented from Graphics Gems IV */
143 
144  Matrix<Scalar, 3, 1> res;
145 
146  const Index odd = ((a0 + 1) % 3 == a1) ? 0 : 1;
147  const Index i = a0;
148  const Index j = (a0 + 1 + odd) % 3;
149  const Index k = (a0 + 2 - odd) % 3;
150 
151  if (a0 == a2)
152  {
153  res[0] = numext::atan2(coeff(j, i), coeff(k, i));
154  if ((odd && res[0] < Scalar(0)) || ((!odd) && res[0] > Scalar(0)))
155  {
156  if (res[0] > Scalar(0))
157  {
158  res[0] -= Scalar(EIGEN_PI);
159  }
160  else
161  {
162  res[0] += Scalar(EIGEN_PI);
163  }
164 
165  Scalar s2 = numext::hypot(coeff(j, i), coeff(k, i));
166  res[1] = -numext::atan2(s2, coeff(i, i));
167  }
168  else
169  {
170  Scalar s2 = numext::hypot(coeff(j, i), coeff(k, i));
171  res[1] = numext::atan2(s2, coeff(i, i));
172  }
173 
174  // With a=(0,1,0), we have i=0; j=1; k=2, and after computing the first two angles,
175  // we can compute their respective rotation, and apply its inverse to M. Since the result must
176  // be a rotation around x, we have:
177  //
178  // c2 s1.s2 c1.s2 1 0 0
179  // 0 c1 -s1 * M = 0 c3 s3
180  // -s2 s1.c2 c1.c2 0 -s3 c3
181  //
182  // Thus: m11.c1 - m21.s1 = c3 & m12.c1 - m22.s1 = s3
183 
184  Scalar s1 = numext::sin(res[0]);
185  Scalar c1 = numext::cos(res[0]);
186  res[2] = numext::atan2(c1 * coeff(j, k) - s1 * coeff(k, k), c1 * coeff(j, j) - s1 * coeff(k, j));
187  }
188  else
189  {
190  res[0] = numext::atan2(coeff(j, k), coeff(k, k));
191  Scalar c2 = numext::hypot(coeff(i, i), coeff(i, j));
192  if ((odd && res[0] < Scalar(0)) || ((!odd) && res[0] > Scalar(0)))
193  {
194  if (res[0] > Scalar(0))
195  {
196  res[0] -= Scalar(EIGEN_PI);
197  }
198  else
199  {
200  res[0] += Scalar(EIGEN_PI);
201  }
202  res[1] = numext::atan2(-coeff(i, k), -c2);
203  }
204  else
205  {
206  res[1] = numext::atan2(-coeff(i, k), c2);
207  }
208  Scalar s1 = numext::sin(res[0]);
209  Scalar c1 = numext::cos(res[0]);
210  res[2] = numext::atan2(s1 * coeff(k, i) - c1 * coeff(j, i), c1 * coeff(j, j) - s1 * coeff(k, j));
211  }
212  if (!odd)
213  {
214  res = -res;
215  }
216 
217  return res;
218 }
#define EIGEN_PI
Definition: MathFunctions.h:16

◆ hnormalized() [1/2]

template<typename Derived >
const MatrixBase< Derived >::HNormalizedReturnType Eigen::MatrixBase< Derived >::hnormalized
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);
FixedSegmentReturnType<... >::Type head(NType n)
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:69

Example:

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;
Array< int, Dynamic, 1 > v
Projective3d P(Matrix4d::Random())
TransposeReturnType transpose()
Definition: Transpose.h:184
static const RandomReturnType Random()
Definition: Random.h:114
Transform< double, 3, Projective > Projective3d
Definition: Transform.h:732
const HNormalizedReturnType hnormalized() const
homogeneous normalization
Definition: Homogeneous.h:176
Matrix< double, 4, 1 > Vector4d
4×1 vector of type double.
Definition: Matrix.h:502

Output:

v                   =   0.68 -0.211  0.566  0.597]^T
v.hnormalized()     =   1.14 -0.354  0.949]^T
P*v                 = 0.663 -0.16 -0.13  0.91]^T
(P*v).hnormalized() =  0.729 -0.176 -0.143]^T
See also
VectorwiseOp::hnormalized()

Definition at line 176 of file Homogeneous.h.

177 {
179  return ConstStartMinusOne(derived(),0,0,
180  ColsAtCompileTime==1?size()-1:1,
181  ColsAtCompileTime==1?1:size()-1) / coeff(size()-1);
182 }
#define EIGEN_STATIC_ASSERT_VECTOR_ONLY(TYPE)
Definition: StaticAssert.h:36
Block< const Derived, internal::traits< Derived >::ColsAtCompileTime==1 ? SizeMinusOne :1, internal::traits< Derived >::ColsAtCompileTime==1 ? 1 :SizeMinusOne > ConstStartMinusOne
Definition: MatrixBase.h:420

◆ hnormalized() [2/2]

template<typename ExpressionType , int Direction>
const VectorwiseOp< ExpressionType, Direction >::HNormalizedReturnType Eigen::VectorwiseOp< ExpressionType, Direction >::hnormalized
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:

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;
Matrix4Xd M
Matrix< double, 4, Dynamic > Matrix4Xd
4×Dynamic matrix of type double.
Definition: Matrix.h:502

Output:

The matrix M is:
   0.68   0.823  -0.444   -0.27   0.271
 -0.211  -0.605   0.108  0.0268   0.435
  0.566   -0.33 -0.0452   0.904  -0.717
  0.597   0.536   0.258   0.832   0.214

M.colwise().hnormalized():
  1.14   1.53  -1.72 -0.325   1.27
-0.354  -1.13  0.419 0.0322   2.03
 0.949 -0.614 -0.175   1.09  -3.35

P*M:
  0.186  -0.589   0.369    1.33   -1.23
 -0.871  -0.337   0.127  -0.715   0.091
 -0.158 -0.0104   0.312   0.429  -0.478
  0.992   0.777  -0.373   0.468  -0.651

(P*M).colwise().hnormalized():
  0.188  -0.759  -0.989    2.85    1.89
 -0.877  -0.433  -0.342   -1.53   -0.14
  -0.16 -0.0134  -0.837   0.915   0.735

See also
MatrixBase::hnormalized()

Definition at line 200 of file Homogeneous.h.

201 {
202  return HNormalized_Block(_expression(),0,0,
203  Direction==Vertical ? _expression().rows()-1 : _expression().rows(),
204  Direction==Horizontal ? _expression().cols()-1 : _expression().cols()).cwiseQuotient(
205  Replicate<HNormalized_Factors,
206  Direction==Vertical ? HNormalized_SizeMinusOne : 1,
207  Direction==Horizontal ? HNormalized_SizeMinusOne : 1>
209  Direction==Vertical ? _expression().rows()-1:0,
210  Direction==Horizontal ? _expression().cols()-1:0,
211  Direction==Vertical ? 1 : _expression().rows(),
212  Direction==Horizontal ? 1 : _expression().cols()),
213  Direction==Vertical ? _expression().rows()-1 : 1,
214  Direction==Horizontal ? _expression().cols()-1 : 1));
215 }
Block< const ExpressionType, Direction==Vertical ? int(HNormalized_SizeMinusOne) :int(internal::traits< ExpressionType >::RowsAtCompileTime), Direction==Horizontal ? int(HNormalized_SizeMinusOne) :int(internal::traits< ExpressionType >::ColsAtCompileTime)> HNormalized_Block
Definition: VectorwiseOp.h:728
Block< const ExpressionType, Direction==Vertical ? 1 :int(internal::traits< ExpressionType >::RowsAtCompileTime), Direction==Horizontal ? 1 :int(internal::traits< ExpressionType >::ColsAtCompileTime)> HNormalized_Factors
Definition: VectorwiseOp.h:732
@ Horizontal
Definition: Constants.h:269

◆ homogeneous() [1/2]

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

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;
EvalReturnType eval() const
Definition: DenseBase.h:405
Matrix< double, 3, 1 > Vector3d
3×1 vector of type double.
Definition: Matrix.h:502

Output:

v                                   = [  0.68 -0.211  0.566]^T
h.homogeneous()                     = [  0.68 -0.211  0.566      1]^T
(P * v.homogeneous())               = [  1.27  0.772 0.0154 -0.419]^T
(P * v.homogeneous()).hnormalized() = [  -3.03   -1.84 -0.0367]^T
See also
VectorwiseOp::homogeneous(), class Homogeneous

Definition at line 134 of file Homogeneous.h.

135 {
137  return HomogeneousReturnType(derived());
138 }
Homogeneous< Derived, HomogeneousReturnTypeDirection > HomogeneousReturnType
Definition: MatrixBase.h:411

◆ homogeneous() [2/2]

template<typename ExpressionType , int Direction>
Homogeneous< ExpressionType, Direction > Eigen::VectorwiseOp< ExpressionType, Direction >::homogeneous
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:

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;
Matrix< double, 3, Dynamic > Matrix3Xd
3×Dynamic matrix of type double.
Definition: Matrix.h:502

Output:

The matrix M is:
   0.68   0.597   -0.33   0.108   -0.27
 -0.211   0.823   0.536 -0.0452  0.0268
  0.566  -0.605  -0.444   0.258   0.904

M.colwise().homogeneous():
   0.68   0.597   -0.33   0.108   -0.27
 -0.211   0.823   0.536 -0.0452  0.0268
  0.566  -0.605  -0.444   0.258   0.904
      1       1       1       1       1

P * M.colwise().homogeneous():
0.0832 -0.477  -1.21 -0.545 -0.452
 0.998  0.779  0.695  0.894  0.277
-0.271 -0.608 -0.895 -0.544 -0.874
-0.728 -0.551  0.202  -0.21 -0.469

P * M.colwise().homogeneous().hnormalized(): 
-0.114  0.866     -6    2.6  0.962
 -1.37  -1.41   3.44  -4.27 -0.591
 0.373    1.1  -4.43    2.6   1.86

See also
MatrixBase::homogeneous(), class Homogeneous

Definition at line 152 of file Homogeneous.h.

153 {
155 }
Homogeneous< ExpressionType, Direction > HomogeneousReturnType
Definition: VectorwiseOp.h:709

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

Definition at line 97 of file Umeyama.h.

98 {
99  typedef typename internal::umeyama_transform_matrix_type<Derived, OtherDerived>::type TransformationMatrixType;
100  typedef typename internal::traits<TransformationMatrixType>::Scalar Scalar;
101  typedef typename NumTraits<Scalar>::Real RealScalar;
102 
103  EIGEN_STATIC_ASSERT(!NumTraits<Scalar>::IsComplex, NUMERIC_TYPE_MUST_BE_REAL)
104  EIGEN_STATIC_ASSERT((internal::is_same<Scalar, typename internal::traits<OtherDerived>::Scalar>::value),
105  YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
106 
107  enum { Dimension = internal::min_size_prefer_dynamic(Derived::RowsAtCompileTime, OtherDerived::RowsAtCompileTime) };
108 
109  typedef Matrix<Scalar, Dimension, 1> VectorType;
110  typedef Matrix<Scalar, Dimension, Dimension> MatrixType;
111  typedef typename internal::plain_matrix_type_row_major<Derived>::type RowMajorMatrixType;
112 
113  const Index m = src.rows(); // dimension
114  const Index n = src.cols(); // number of measurements
115 
116  // required for demeaning ...
117  const RealScalar one_over_n = RealScalar(1) / static_cast<RealScalar>(n);
118 
119  // computation of mean
120  const VectorType src_mean = src.rowwise().sum() * one_over_n;
121  const VectorType dst_mean = dst.rowwise().sum() * one_over_n;
122 
123  // demeaning of src and dst points
124  const RowMajorMatrixType src_demean = src.colwise() - src_mean;
125  const RowMajorMatrixType dst_demean = dst.colwise() - dst_mean;
126 
127  // Eq. (38)
128  const MatrixType sigma = one_over_n * dst_demean * src_demean.transpose();
129 
130  JacobiSVD<MatrixType, ComputeFullU | ComputeFullV> svd(sigma);
131 
132  // Initialize the resulting transformation with an identity matrix...
133  TransformationMatrixType Rt = TransformationMatrixType::Identity(m+1,m+1);
134 
135  // Eq. (39)
136  VectorType S = VectorType::Ones(m);
137 
138  if ( svd.matrixU().determinant() * svd.matrixV().determinant() < 0 ) {
139  Index tmp = m - 1;
140  S(tmp) = -1;
141  }
142 
143  // Eq. (40) and (43)
144  Rt.block(0,0,m,m).noalias() = svd.matrixU() * S.asDiagonal() * svd.matrixV().transpose();
145 
146  if (with_scaling)
147  {
148  // Eq. (36)-(37)
149  const Scalar src_var = src_demean.rowwise().squaredNorm().sum() * one_over_n;
150 
151  // Eq. (42)
152  const Scalar c = Scalar(1)/src_var * svd.singularValues().dot(S);
153 
154  // Eq. (41)
155  Rt.col(m).head(m) = dst_mean;
156  Rt.col(m).head(m).noalias() -= c*Rt.topLeftCorner(m,m)*src_mean;
157  Rt.block(0,0,m,m) *= c;
158  }
159  else
160  {
161  Rt.col(m).head(m) = dst_mean;
162  Rt.col(m).head(m).noalias() -= Rt.topLeftCorner(m,m)*src_mean;
163  }
164 
165  return Rt;
166 }
Matrix3f m
int n
Array33i c
cout<< "Here is the matrix m:"<< endl<< m<< endl;JacobiSVD< MatrixXf, ComputeThinU|ComputeThinV > svd(m)
Matrix< float, 1, Dynamic > MatrixType
constexpr int min_size_prefer_dynamic(A a, B b)
Definition: Meta.h:537
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82

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

Definition at line 273 of file OrthoMethods.h.

274 {
276  return internal::unitOrthogonal_selector<Derived>::run(derived());
277 }