mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-01-06 14:14:46 +08:00
231 lines
8.3 KiB
Plaintext
231 lines
8.3 KiB
Plaintext
namespace Eigen {
|
|
|
|
/** \page TutorialGeometry Tutorial 2/3 - Geometry
|
|
\ingroup Tutorial
|
|
|
|
<div class="eimainmenu">\ref index "Overview"
|
|
| \ref TutorialCore "Core features"
|
|
| \b Geometry
|
|
| \ref TutorialAdvancedLinearAlgebra "Advanced linear algebra"
|
|
</div>
|
|
|
|
In this tutorial chapter we will shortly introduce the many possibilities offered by the \ref GeometryModule "geometry module",
|
|
namely 2D and 3D rotations and affine transformations.
|
|
|
|
\b Table \b of \b contents
|
|
- \ref TutorialGeoElementaryTransformations
|
|
- \ref TutorialGeoCommontransformationAPI
|
|
- \ref TutorialGeoTransform
|
|
- \ref TutorialGeoEulerAngles
|
|
|
|
\section TutorialGeoElementaryTransformations Transformation types
|
|
|
|
<table class="tutorial_code">
|
|
<tr><td>Transformation type</td><td>Typical initialization code</td></tr>
|
|
<tr><td>
|
|
\ref Rotation2D "2D rotation" from an angle</td><td>\code
|
|
Rotation2D<float> rot2(angle_in_radian);\endcode</td></tr>
|
|
<tr><td>
|
|
3D rotation as an \ref AngleAxis "angle + axis"</td><td>\code
|
|
AngleAxis<float> aa(angle_in_radian, Vector3f(ax,ay,az));\endcode</td></tr>
|
|
<tr><td>
|
|
3D rotation as a \ref Quaternion "quaternion"</td><td>\code
|
|
Quaternion<float> q = AngleAxis<float>(angle_in_radian, axis);\endcode</td></tr>
|
|
<tr><td>
|
|
N-D Scaling</td><td>\code
|
|
Scaling<float,2>(sx, sy)
|
|
Scaling<float,3>(sx, sy, sz)
|
|
Scaling<float,N>(s)
|
|
Scaling<float,N>(vecN)\endcode</td></tr>
|
|
<tr><td>
|
|
N-D Translation</td><td>\code
|
|
Translation<float,2>(tx, ty)
|
|
Translation<float,3>(tx, ty, tz)
|
|
Translation<float,N>(s)
|
|
Translation<float,N>(vecN)\endcode</td></tr>
|
|
<tr><td>
|
|
N-D \ref TutorialGeoTransform "Affine transformation"</td><td>\code
|
|
Transform<float,N> t = concatenation_of_any_transformations;
|
|
Transform<float,3> t = Translation3f(p) * AngleAxisf(a,axis) * Scaling3f(s);\endcode</td></tr>
|
|
<tr><td>
|
|
N-D Linear transformations \n
|
|
<em class=note>(pure rotations, \n scaling, etc.)</em></td><td>\code
|
|
Matrix<float,N> t = concatenation_of_rotations_and_scalings;
|
|
Matrix<float,2> t = Rotation2Df(a) * Scaling2f(s);
|
|
Matrix<float,3> t = AngleAxisf(a,axis) * Scaling3f(s);\endcode</td></tr>
|
|
</table>
|
|
|
|
<strong>Notes on rotations</strong>\n To transform more than a single vector the preferred
|
|
representations are rotation matrices, while for other usages Quaternion is the
|
|
representation of choice as they are compact, fast and stable. Finally Rotation2D and
|
|
AngleAxis are mainly convenient types to create other rotation objects.
|
|
|
|
<strong>Notes on Translation and Scaling</strong>\n Likewise AngleAxis, these classes were
|
|
designed to simplify the creation/initialization of linear (Matrix) and affine (Transform)
|
|
transformations. Nevertheless, unlike AngleAxis which is inefficient to use, these classes
|
|
might still be interesting to write generic and efficient algorithms taking as input any
|
|
kind of transformations.
|
|
|
|
Any of the above transformation types can be converted to any other types of the same nature,
|
|
or to a more generic type. Here are come additional examples:
|
|
<table class="tutorial_code">
|
|
<tr><td>\code
|
|
Rotation2Df r = Matrix2f(..); // assumes a pure rotation matrix
|
|
AngleAxisf aa = Quaternionf(..);
|
|
AngleAxisf aa = Matrix3f(..); // assumes a pure rotation matrix
|
|
Matrix2f m = Rotation2Df(..);
|
|
Matrix3f m = Quaternionf(..); Matrix3f m = Scaling3f(..);
|
|
Transform3f m = AngleAxis3f(..); Transform3f m = Scaling3f(..);
|
|
Transform3f m = Translation3f(..); Transform3f m = Matrix3f(..);
|
|
\endcode</td></tr>
|
|
</table>
|
|
|
|
|
|
<a href="#" class="top">top</a>\section TutorialGeoCommontransformationAPI Common API across transformation types
|
|
|
|
To some extent, Eigen's \ref GeometryModule "geometry module" allows you to write
|
|
generic algorithms working on any kind of transformation representations:
|
|
<table class="tutorial_code">
|
|
<tr><td>
|
|
Concatenation of two transformations</td><td>\code
|
|
gen1 * gen2;\endcode</td></tr>
|
|
<tr><td>Apply the transformation to a vector</td><td>\code
|
|
vec2 = gen1 * vec1;\endcode</td></tr>
|
|
<tr><td>Get the inverse of the transformation</td><td>\code
|
|
gen2 = gen1.inverse();\endcode</td></tr>
|
|
<tr><td>Spherical interpolation \n (Rotation2D and Quaternion only)</td><td>\code
|
|
rot3 = rot1.slerp(alpha,rot2);\endcode</td></tr>
|
|
</table>
|
|
|
|
|
|
|
|
<a href="#" class="top">top</a>\section TutorialGeoTransform Affine transformations
|
|
Generic affine transformations are represented by the Transform class which internaly
|
|
is a (Dim+1)^2 matrix. In Eigen we have chosen to not distinghish between points and
|
|
vectors such that all points are actually represented by displacement vectors from the
|
|
origin ( \f$ \mathbf{p} \equiv \mathbf{p}-0 \f$ ). With that in mind, real points and
|
|
vector distinguish when the transformation is applied.
|
|
<table class="tutorial_code">
|
|
<tr><td>
|
|
Apply the transformation to a \b point </td><td>\code
|
|
VectorNf p1, p2;
|
|
p2 = t * p1;\endcode</td></tr>
|
|
<tr><td>
|
|
Apply the transformation to a \b vector </td><td>\code
|
|
VectorNf vec1, vec2;
|
|
vec2 = t.linear() * vec1;\endcode</td></tr>
|
|
<tr><td>
|
|
Apply a \em general transformation \n to a \b normal \b vector
|
|
(<a href="http://www.cgafaq.info/wiki/Transforming_normals">explanations</a>)</td><td>\code
|
|
VectorNf n1, n2;
|
|
MatrixNf normalMatrix = t.linear().inverse().transpose();
|
|
n2 = (normalMatrix * n1).normalized();\endcode</td></tr>
|
|
<tr><td>
|
|
Apply a transformation with \em pure \em rotation \n to a \b normal \b vector
|
|
(no scaling, no shear)</td><td>\code
|
|
n2 = t.linear() * n1;\endcode</td></tr>
|
|
<tr><td>
|
|
OpenGL compatibility \b 3D </td><td>\code
|
|
glLoadMatrixf(t.data());\endcode</td></tr>
|
|
<tr><td>
|
|
OpenGL compatibility \b 2D </td><td>\code
|
|
Transform3f aux(Transform3f::Identity);
|
|
aux.linear().corner<2,2>(TopLeft) = t.linear();
|
|
aux.translation().start<2>() = t.translation();
|
|
glLoadMatrixf(aux.data());\endcode</td></tr>
|
|
</table>
|
|
|
|
\b Component \b accessors</td></tr>
|
|
<table class="tutorial_code">
|
|
<tr><td>
|
|
full read-write access to the internal matrix</td><td>\code
|
|
t.matrix() = matN1xN1; // N1 means N+1
|
|
matN1xN1 = t.matrix();
|
|
\endcode</td></tr>
|
|
<tr><td>
|
|
coefficient accessors</td><td>\code
|
|
t(i,j) = scalar; <=> t.matrix()(i,j) = scalar;
|
|
scalar = t(i,j); <=> scalar = t.matrix()(i,j);
|
|
\endcode</td></tr>
|
|
<tr><td>
|
|
translation part</td><td>\code
|
|
t.translation() = vecN;
|
|
vecN = t.translation();
|
|
\endcode</td></tr>
|
|
<tr><td>
|
|
linear part</td><td>\code
|
|
t.linear() = matNxN;
|
|
matNxN = t.linear();
|
|
\endcode</td></tr>
|
|
<tr><td>
|
|
extract the rotation matrix</td><td>\code
|
|
matNxN = t.extractRotation();
|
|
\endcode</td></tr>
|
|
</table>
|
|
|
|
|
|
\b Transformation \b creation \n
|
|
While transformation objects can be created and updated concatenating elementary transformations,
|
|
the Transform class also features a procedural API:
|
|
<table class="tutorial_code">
|
|
<tr><td></td><td>\b procedurale \b API </td><td>\b equivalent \b natural \b API </td></tr>
|
|
<tr><td>Translation</td><td>\code
|
|
t.translate(Vector_(tx,ty,..));
|
|
t.pretranslate(Vector_(tx,ty,..));
|
|
\endcode</td><td>\code
|
|
t *= Translation_(tx,ty,..);
|
|
t = Translation_(tx,ty,..) * t;
|
|
\endcode</td></tr>
|
|
<tr><td>\b Rotation \n <em class="note">In 2D, any_rotation can also \n be an angle in radian</em></td><td>\code
|
|
t.rotate(any_rotation);
|
|
t.prerotate(any_rotation);
|
|
\endcode</td><td>\code
|
|
t *= any_rotation;
|
|
t = any_rotation * t;
|
|
\endcode</td></tr>
|
|
<tr><td>Scaling</td><td>\code
|
|
t.scale(Vector_(sx,sy,..));
|
|
t.scale(s);
|
|
t.prescale(Vector_(sx,sy,..));
|
|
t.prescale(s);
|
|
\endcode</td><td>\code
|
|
t *= Scaling_(sx,sy,..);
|
|
t *= Scaling_(s);
|
|
t = Scaling_(sx,sy,..) * t;
|
|
t = Scaling_(s) * t;
|
|
\endcode</td></tr>
|
|
<tr><td>Shear transformation \n ( \b 2D \b only ! )</td><td>\code
|
|
t.shear(sx,sy);
|
|
t.preshear(sx,sy);
|
|
\endcode</td><td></td></tr>
|
|
</table>
|
|
|
|
Note that in both API, any many transformations can be concatenated in a single expression as shown in the two following equivalent examples:
|
|
<table class="tutorial_code">
|
|
<tr><td>\code
|
|
t.pretranslate(..).rotate(..).translate(..).scale(..);
|
|
\endcode</td></tr>
|
|
<tr><td>\code
|
|
t = Translation_(..) * t * RotationType(..) * Translation_(..) * Scaling_(..);
|
|
\endcode</td></tr>
|
|
</table>
|
|
|
|
|
|
|
|
<a href="#" class="top">top</a>\section TutorialGeoEulerAngles Euler angles
|
|
<table class="tutorial_code">
|
|
<tr><td style="max-width:30em;">
|
|
Euler angles might be convenient to create rotation objects.
|
|
On the other hand, since there exist 24 differents convensions,they are pretty confusing to use. This example shows how
|
|
to create a rotation matrix according to the 2-1-2 convention.</td><td>\code
|
|
Matrix3f m;
|
|
m = AngleAxisf(angle1, Vector3f::UnitZ())
|
|
* * AngleAxisf(angle2, Vector3f::UnitY())
|
|
* * AngleAxisf(angle3, Vector3f::UnitZ());
|
|
\endcode</td></tr>
|
|
</table>
|
|
|
|
*/
|
|
|
|
}
|