1 namespace Eigen { 2 3 /** \eigenManualPage TutorialGeometry Space transformations 4 5 In this page, we will introduce the many possibilities offered by the \ref Geometry_Module "geometry module" to deal with 2D and 3D rotations and projective or affine transformations. 6 7 \eigenAutoToc 8 9 Eigen's Geometry module provides two different kinds of geometric transformations: 10 - Abstract transformations, such as rotations (represented by \ref AngleAxis "angle and axis" or by a \ref Quaternion "quaternion"), \ref Translation "translations", \ref Scaling "scalings". These transformations are NOT represented as matrices, but you can nevertheless mix them with matrices and vectors in expressions, and convert them to matrices if you wish. 11 - Projective or affine transformation matrices: see the Transform class. These are really matrices. 12 13 \note If you are working with OpenGL 4x4 matrices then Affine3f and Affine3d are what you want. Since Eigen defaults to column-major storage, you can directly use the Transform::data() method to pass your transformation matrix to OpenGL. 14 15 You can construct a Transform from an abstract transformation, like this: 16 \code 17 Transform t(AngleAxis(angle,axis)); 18 \endcode 19 or like this: 20 \code 21 Transform t; 22 t = AngleAxis(angle,axis); 23 \endcode 24 But note that unfortunately, because of how C++ works, you can \b not do this: 25 \code 26 Transform t = AngleAxis(angle,axis); 27 \endcode 28 <span class="note">\b Explanation: In the C++ language, this would require Transform to have a non-explicit conversion constructor from AngleAxis, but we really don't want to allow implicit casting here. 29 </span> 30 31 \section TutorialGeoElementaryTransformations Transformation types 32 33 <table class="manual"> 34 <tr><th>Transformation type</th><th>Typical initialization code</th></tr> 35 <tr><td> 36 \ref Rotation2D "2D rotation" from an angle</td><td>\code 37 Rotation2D<float> rot2(angle_in_radian);\endcode</td></tr> 38 <tr class="alt"><td> 39 3D rotation as an \ref AngleAxis "angle + axis"</td><td>\code 40 AngleAxis<float> aa(angle_in_radian, Vector3f(ax,ay,az));\endcode 41 <span class="note">The axis vector must be normalized.</span></td></tr> 42 <tr><td> 43 3D rotation as a \ref Quaternion "quaternion"</td><td>\code 44 Quaternion<float> q; q = AngleAxis<float>(angle_in_radian, axis);\endcode</td></tr> 45 <tr class="alt"><td> 46 N-D Scaling</td><td>\code 47 Scaling(sx, sy) 48 Scaling(sx, sy, sz) 49 Scaling(s) 50 Scaling(vecN)\endcode</td></tr> 51 <tr><td> 52 N-D Translation</td><td>\code 53 Translation<float,2>(tx, ty) 54 Translation<float,3>(tx, ty, tz) 55 Translation<float,N>(s) 56 Translation<float,N>(vecN)\endcode</td></tr> 57 <tr class="alt"><td> 58 N-D \ref TutorialGeoTransform "Affine transformation"</td><td>\code 59 Transform<float,N,Affine> t = concatenation_of_any_transformations; 60 Transform<float,3,Affine> t = Translation3f(p) * AngleAxisf(a,axis) * Scaling(s);\endcode</td></tr> 61 <tr><td> 62 N-D Linear transformations \n 63 <em class=note>(pure rotations, \n scaling, etc.)</em></td><td>\code 64 Matrix<float,N> t = concatenation_of_rotations_and_scalings; 65 Matrix<float,2> t = Rotation2Df(a) * Scaling(s); 66 Matrix<float,3> t = AngleAxisf(a,axis) * Scaling(s);\endcode</td></tr> 67 </table> 68 69 <strong>Notes on rotations</strong>\n To transform more than a single vector the preferred 70 representations are rotation matrices, while for other usages Quaternion is the 71 representation of choice as they are compact, fast and stable. Finally Rotation2D and 72 AngleAxis are mainly convenient types to create other rotation objects. 73 74 <strong>Notes on Translation and Scaling</strong>\n Like AngleAxis, these classes were 75 designed to simplify the creation/initialization of linear (Matrix) and affine (Transform) 76 transformations. Nevertheless, unlike AngleAxis which is inefficient to use, these classes 77 might still be interesting to write generic and efficient algorithms taking as input any 78 kind of transformations. 79 80 Any of the above transformation types can be converted to any other types of the same nature, 81 or to a more generic type. Here are some additional examples: 82 <table class="manual"> 83 <tr><td>\code 84 Rotation2Df r; r = Matrix2f(..); // assumes a pure rotation matrix 85 AngleAxisf aa; aa = Quaternionf(..); 86 AngleAxisf aa; aa = Matrix3f(..); // assumes a pure rotation matrix 87 Matrix2f m; m = Rotation2Df(..); 88 Matrix3f m; m = Quaternionf(..); Matrix3f m; m = Scaling(..); 89 Affine3f m; m = AngleAxis3f(..); Affine3f m; m = Scaling(..); 90 Affine3f m; m = Translation3f(..); Affine3f m; m = Matrix3f(..); 91 \endcode</td></tr> 92 </table> 93 94 95 <a href="#" class="top">top</a>\section TutorialGeoCommontransformationAPI Common API across transformation types 96 97 To some extent, Eigen's \ref Geometry_Module "geometry module" allows you to write 98 generic algorithms working on any kind of transformation representations: 99 <table class="manual"> 100 <tr><td> 101 Concatenation of two transformations</td><td>\code 102 gen1 * gen2;\endcode</td></tr> 103 <tr class="alt"><td>Apply the transformation to a vector</td><td>\code 104 vec2 = gen1 * vec1;\endcode</td></tr> 105 <tr><td>Get the inverse of the transformation</td><td>\code 106 gen2 = gen1.inverse();\endcode</td></tr> 107 <tr class="alt"><td>Spherical interpolation \n (Rotation2D and Quaternion only)</td><td>\code 108 rot3 = rot1.slerp(alpha,rot2);\endcode</td></tr> 109 </table> 110 111 112 113 <a href="#" class="top">top</a>\section TutorialGeoTransform Affine transformations 114 Generic affine transformations are represented by the Transform class which internaly 115 is a (Dim+1)^2 matrix. In Eigen we have chosen to not distinghish between points and 116 vectors such that all points are actually represented by displacement vectors from the 117 origin ( \f$ \mathbf{p} \equiv \mathbf{p}-0 \f$ ). With that in mind, real points and 118 vector distinguish when the transformation is applied. 119 <table class="manual"> 120 <tr><td> 121 Apply the transformation to a \b point </td><td>\code 122 VectorNf p1, p2; 123 p2 = t * p1;\endcode</td></tr> 124 <tr class="alt"><td> 125 Apply the transformation to a \b vector </td><td>\code 126 VectorNf vec1, vec2; 127 vec2 = t.linear() * vec1;\endcode</td></tr> 128 <tr><td> 129 Apply a \em general transformation \n to a \b normal \b vector \n 130 </td><td>\code 131 VectorNf n1, n2; 132 MatrixNf normalMatrix = t.linear().inverse().transpose(); 133 n2 = (normalMatrix * n1).normalized();\endcode</td></tr> 134 <tr><td colspan="2">(See subject 5.27 of this <a href="http://www.faqs.org/faqs/graphics/algorithms-faq">faq</a> for the explanations)</td></tr> 135 <tr class="alt"><td> 136 Apply a transformation with \em pure \em rotation \n to a \b normal \b vector 137 (no scaling, no shear)</td><td>\code 138 n2 = t.linear() * n1;\endcode</td></tr> 139 <tr><td> 140 OpenGL compatibility \b 3D </td><td>\code 141 glLoadMatrixf(t.data());\endcode</td></tr> 142 <tr class="alt"><td> 143 OpenGL compatibility \b 2D </td><td>\code 144 Affine3f aux(Affine3f::Identity()); 145 aux.linear().topLeftCorner<2,2>() = t.linear(); 146 aux.translation().start<2>() = t.translation(); 147 glLoadMatrixf(aux.data());\endcode</td></tr> 148 </table> 149 150 \b Component \b accessors 151 <table class="manual"> 152 <tr><td> 153 full read-write access to the internal matrix</td><td>\code 154 t.matrix() = matN1xN1; // N1 means N+1 155 matN1xN1 = t.matrix(); 156 \endcode</td></tr> 157 <tr class="alt"><td> 158 coefficient accessors</td><td>\code 159 t(i,j) = scalar; <=> t.matrix()(i,j) = scalar; 160 scalar = t(i,j); <=> scalar = t.matrix()(i,j); 161 \endcode</td></tr> 162 <tr><td> 163 translation part</td><td>\code 164 t.translation() = vecN; 165 vecN = t.translation(); 166 \endcode</td></tr> 167 <tr class="alt"><td> 168 linear part</td><td>\code 169 t.linear() = matNxN; 170 matNxN = t.linear(); 171 \endcode</td></tr> 172 <tr><td> 173 extract the rotation matrix</td><td>\code 174 matNxN = t.rotation(); 175 \endcode</td></tr> 176 </table> 177 178 179 \b Transformation \b creation \n 180 While transformation objects can be created and updated concatenating elementary transformations, 181 the Transform class also features a procedural API: 182 <table class="manual"> 183 <tr><th></th><th>procedural API</th><th>equivalent natural API </th></tr> 184 <tr><td>Translation</td><td>\code 185 t.translate(Vector_(tx,ty,..)); 186 t.pretranslate(Vector_(tx,ty,..)); 187 \endcode</td><td>\code 188 t *= Translation_(tx,ty,..); 189 t = Translation_(tx,ty,..) * t; 190 \endcode</td></tr> 191 <tr class="alt"><td>\b Rotation \n <em class="note">In 2D and for the procedural API, any_rotation can also \n be an angle in radian</em></td><td>\code 192 t.rotate(any_rotation); 193 t.prerotate(any_rotation); 194 \endcode</td><td>\code 195 t *= any_rotation; 196 t = any_rotation * t; 197 \endcode</td></tr> 198 <tr><td>Scaling</td><td>\code 199 t.scale(Vector_(sx,sy,..)); 200 t.scale(s); 201 t.prescale(Vector_(sx,sy,..)); 202 t.prescale(s); 203 \endcode</td><td>\code 204 t *= Scaling(sx,sy,..); 205 t *= Scaling(s); 206 t = Scaling(sx,sy,..) * t; 207 t = Scaling(s) * t; 208 \endcode</td></tr> 209 <tr class="alt"><td>Shear transformation \n ( \b 2D \b only ! )</td><td>\code 210 t.shear(sx,sy); 211 t.preshear(sx,sy); 212 \endcode</td><td></td></tr> 213 </table> 214 215 Note that in both API, any many transformations can be concatenated in a single expression as shown in the two following equivalent examples: 216 <table class="manual"> 217 <tr><td>\code 218 t.pretranslate(..).rotate(..).translate(..).scale(..); 219 \endcode</td></tr> 220 <tr><td>\code 221 t = Translation_(..) * t * RotationType(..) * Translation_(..) * Scaling(..); 222 \endcode</td></tr> 223 </table> 224 225 226 227 <a href="#" class="top">top</a>\section TutorialGeoEulerAngles Euler angles 228 <table class="manual"> 229 <tr><td style="max-width:30em;"> 230 Euler angles might be convenient to create rotation objects. 231 On the other hand, since there exist 24 different conventions, they are pretty confusing to use. This example shows how 232 to create a rotation matrix according to the 2-1-2 convention.</td><td>\code 233 Matrix3f m; 234 m = AngleAxisf(angle1, Vector3f::UnitZ()) 235 * AngleAxisf(angle2, Vector3f::UnitY()) 236 * AngleAxisf(angle3, Vector3f::UnitZ()); 237 \endcode</td></tr> 238 </table> 239 240 */ 241 242 } 243