aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/src/Geometry/Hyperplane.h22
-rw-r--r--Eigen/src/Geometry/Transform.h46
-rw-r--r--doc/QuickStartGuide.dox221
-rw-r--r--doc/TutorialGeometry.dox6
-rw-r--r--doc/eigendoxy.css37
5 files changed, 84 insertions, 248 deletions
diff --git a/Eigen/src/Geometry/Hyperplane.h b/Eigen/src/Geometry/Hyperplane.h
index 5a8657f31..936808ec2 100644
--- a/Eigen/src/Geometry/Hyperplane.h
+++ b/Eigen/src/Geometry/Hyperplane.h
@@ -23,8 +23,8 @@
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.
-#ifndef EIGEN_Hyperplane_H
-#define EIGEN_Hyperplane_H
+#ifndef EIGEN_HYPERPLANE_H
+#define EIGEN_HYPERPLANE_H
/** \geometry_module \ingroup GeometryModule
*
@@ -279,20 +279,4 @@ inline _Scalar ParametrizedLine<_Scalar, _AmbientDim>::intersection(const Hyperp
/(direction().dot(hyperplane.normal()));
}
-#if 0
-/** \addtogroup GeometryModule */
-//@{
-typedef Hyperplane<float, 2> Hyperplane2f;
-typedef Hyperplane<double,2> Hyperplane2d;
-typedef Hyperplane<float, 3> Hyperplane3f;
-typedef Hyperplane<double,3> Hyperplane3d;
-
-typedef Hyperplane<float, 3> Planef;
-typedef Hyperplane<double,3> Planed;
-
-typedef Hyperplane<float, Dynamic> HyperplaneXf;
-typedef Hyperplane<double,Dynamic> HyperplaneXd;
-//@}
-#endif
-
-#endif // EIGEN_Hyperplane_H
+#endif // EIGEN_HYPERPLANE_H
diff --git a/Eigen/src/Geometry/Transform.h b/Eigen/src/Geometry/Transform.h
index b003f20ed..3586fb60f 100644
--- a/Eigen/src/Geometry/Transform.h
+++ b/Eigen/src/Geometry/Transform.h
@@ -55,8 +55,8 @@ struct ei_transform_product_impl;
* The homography is internally represented and stored as a (Dim+1)^2 matrix which
* is available through the matrix() method.
*
- * Conversion methods from/to Qt's QMatrix are available if the preprocessor token
- * EIGEN_QT_SUPPORT is defined.
+ * Conversion methods from/to Qt's QMatrix and QTransform are available if the
+ * preprocessor token EIGEN_QT_SUPPORT is defined.
*
* \sa class Matrix, class Quaternion
*/
@@ -137,6 +137,9 @@ public:
inline Transform(const QMatrix& other);
inline Transform& operator=(const QMatrix& other);
inline QMatrix toQMatrix(void) const;
+ inline Transform(const QTransform& other);
+ inline Transform& operator=(const QTransform& other);
+ inline QTransform toQTransform(void) const;
#endif
/** shortcut for m_matrix(row,col);
@@ -283,6 +286,8 @@ Transform<Scalar,Dim>& Transform<Scalar,Dim>::operator=(const QMatrix& other)
/** \returns a QMatrix from \c *this assuming the dimension is 2.
*
+ * \warning this convertion might loss data if \c *this is not affine
+ *
* This function is available only if the token EIGEN_QT_SUPPORT is defined.
*/
template<typename Scalar, int Dim>
@@ -293,6 +298,43 @@ QMatrix Transform<Scalar,Dim>::toQMatrix(void) const
other.coeffRef(0,1), other.coeffRef(1,1),
other.coeffRef(0,2), other.coeffRef(1,2));
}
+
+/** Initialises \c *this from a QTransform assuming the dimension is 2.
+ *
+ * This function is available only if the token EIGEN_QT_SUPPORT is defined.
+ */
+template<typename Scalar, int Dim>
+Transform<Scalar,Dim>::Transform(const QTransform& other)
+{
+ *this = other;
+}
+
+/** Set \c *this from a QTransform assuming the dimension is 2.
+ *
+ * This function is available only if the token EIGEN_QT_SUPPORT is defined.
+ */
+template<typename Scalar, int Dim>
+Transform<Scalar,Dim>& Transform<Scalar,Dim>::operator=(const QTransform& other)
+{
+ EIGEN_STATIC_ASSERT(Dim==2, you_did_a_programming_error);
+ m_matrix << other.m11(), other.m21(), other.dx(),
+ other.m12(), other.m22(), other.dy(),
+ other.m13(), other.m23(), other.m33();
+ return *this;
+}
+
+/** \returns a QTransform from \c *this assuming the dimension is 2.
+ *
+ * This function is available only if the token EIGEN_QT_SUPPORT is defined.
+ */
+template<typename Scalar, int Dim>
+QMatrix Transform<Scalar,Dim>::toQTransform(void) const
+{
+ EIGEN_STATIC_ASSERT(Dim==2, you_did_a_programming_error);
+ return QTransform(other.coeffRef(0,0), other.coeffRef(1,0), other.coeffRef(2,0)
+ other.coeffRef(0,1), other.coeffRef(1,1), other.coeffRef(2,1)
+ other.coeffRef(0,2), other.coeffRef(1,2), other.coeffRef(2,2);
+}
#endif
/*********************
diff --git a/doc/QuickStartGuide.dox b/doc/QuickStartGuide.dox
index 96497073d..8d1799fdd 100644
--- a/doc/QuickStartGuide.dox
+++ b/doc/QuickStartGuide.dox
@@ -188,7 +188,7 @@ MatrixXi mat2x2 = Map<MatrixXi>(data,2,2);
\subsection TutorialCommaInit Comma initializer
-Eigen also offer a \ref MatrixBaseCommaInitRef "comma initializer syntax" which allows you to set all the coefficients of a matrix to specific values:
+Eigen also offers a \ref MatrixBaseCommaInitRef "comma initializer syntax" which allows you to set all the coefficients of a matrix to specific values:
<table class="tutorial_code"><tr><td>
\include Tutorial_commainit_01.cpp
</td>
@@ -259,7 +259,7 @@ vec3 = vec1.cross(vec2);\endcode</td></tr>
In Eigen only mathematically well defined operators can be used right away,
but don't worry, thanks to the \link Cwise .cwise() \endlink operator prefix,
-Eigen's matrices also provide a very powerful numerical container supporting
+Eigen's matrices are also very powerful as a numerical container supporting
most common coefficient wise operators:
<table class="noborder">
<tr><td>
@@ -322,14 +322,14 @@ mat3 = mat1.cwise().abs2(mat2);
</table>
</td></tr></table>
-<span class="note">\b Side \b note: If you feel the \c .cwise() syntax is too verbose for your taste and don't bother to have non mathematical operator directly available feel free to extend MatrixBase as described \ref ExtendingMatrixBase "here".</span>
+<span class="note">\b Side \b note: If you feel the \c .cwise() syntax is too verbose for your taste and don't bother to have non mathematical operator directly available, then feel free to extend MatrixBase as described \ref ExtendingMatrixBase "here".</span>
<a href="#" class="top">top</a>\section TutorialCoreReductions Reductions
-Eigen provides several several reduction methods such as:
+Eigen provides several reduction methods such as:
\link MatrixBase::minCoeff() minCoeff() \endlink, \link MatrixBase::maxCoeff() maxCoeff() \endlink,
\link MatrixBase::sum() sum() \endlink, \link MatrixBase::trace() trace() \endlink,
\link MatrixBase::norm() norm() \endlink, \link MatrixBase::norm2() norm2() \endlink,
@@ -370,28 +370,22 @@ Read-write access to sub-vectors:
<table class="tutorial_code">
<tr>
<td>Default versions</td>
-<td>Optimized versions when the size is known at compile time</td></tr>
+<td>Optimized versions when the size \n is known at compile time</td></tr>
<td></td>
<tr><td>\code vec1.start(n)\endcode</td><td>\code vec1.start<n>()\endcode</td><td>the first \c n coeffs </td></tr>
<tr><td>\code vec1.end(n)\endcode</td><td>\code vec1.end<n>()\endcode</td><td>the last \c n coeffs </td></tr>
<tr><td>\code vec1.block(pos,n)\endcode</td><td>\code vec1.block<n>(pos)\endcode</td>
- <td>the \c size coeffs in the range [\c pos : \c pos + \c n [</td></tr>
-</table>
-
-Read-write access to sub-matrices:
-<table class="tutorial_code">
-<tr><td>Default versions</td>
-<td>Optimized versions when the size is known at compile time</td><td></td></tr>
+ <td>the \c size coeffs in \n the range [\c pos : \c pos + \c n [</td></tr>
+<tr style="border-style: dashed none dashed none;"><td>
+Read-write access to sub-matrices:</td><td></td><td></td></tr>
<tr>
<td>\code mat1.block(i,j,rows,cols)\endcode
\link MatrixBase::block(int,int,int,int) (more) \endlink</td>
<td>\code mat1.block<rows,cols>(i,j)\endcode
\link MatrixBase::block(int,int) (more) \endlink</td>
- <td>the \c rows x \c cols sub-matrix starting from position (\c i,\c j) </td>
-</tr>
-<tr>
+ <td>the \c rows x \c cols sub-matrix \n starting from position (\c i,\c j)</td></tr><tr>
<td>\code
mat1.corner(TopLeft,rows,cols)
mat1.corner(TopRight,rows,cols)
@@ -493,203 +487,6 @@ forces immediate evaluation of the transpose</td></tr>
-
-
-/** \page TutorialGeometry Tutorial 2/3 - Geometry
- \ingroup Tutorial
- \internal
-
-<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 TutorialGeoRotations
- - \ref TutorialGeoTransformation
-
-\section TutorialGeoRotations 2D and 3D Rotations
-
-\subsection TutorialGeoRotationTypes Rotation types
-
-
-<table class="tutorial_code">
-<tr><td>Rotation type</td><td>Typical initialization code</td><td>Recommended usage</td></tr>
-<tr><td>2D rotation from an angle</td><td>\code
-Rotation2D<float> rot2(angle_in_radian);\endcode</td><td></td></tr>
-<tr><td>2D rotation matrix</td><td>\code
-Matrix2f rotmat2 = Rotation2Df(angle_in_radian);\endcode</td><td></td></tr>
-<tr><td>3D rotation as an angle + axis</td><td>\code
-AngleAxis<float> aa(angle_in_radian, Vector3f(ax,ay,az));\endcode</td><td></td></tr>
-<tr><td>3D rotation as a quaternion</td><td>\code
-Quaternion<float> q = AngleAxis<float>(angle_in_radian, axis);\endcode</td><td></td></tr>
-<tr><td>3D rotation matrix</td><td>\code
-Matrix3f rotmat3 = AngleAxis<float>(angle_in_radian, axis);\endcode</td><td></td></tr>
-</table>
-
-To transform more than a single vector the prefered representations are rotation matrices,
-for other usage Rotation2D and Quaternion are the representations of choice as they are
-more compact, fast and stable. AngleAxis are only useful to create other rotation objects.
-
-\subsection TutorialGeoCommonRotationAPI Common API across rotation types
-
-To some extent, Eigen's \ref GeometryModule "geometry module" allows you to write
-generic algorithms working on both 2D and 3D rotations of any of the five above types.
-The following operation are supported:
-<table class="tutorial_code">
-<tr><td>Convertion from and to any types (of same space dimension)</td><td>\code
-RotType2 a = RotType1();\endcode</td></tr>
-<tr><td>Concatenation of two rotations</td><td>\code
-rot3 = rot1 * rot2;\endcode</td></tr>
-<tr><td>Apply the rotation to a vector</td><td>\code
-vec2 = rot1 * vec1;\endcode</td></tr>
-<tr><td>Get the inverse rotation \n (not always the most effient choice)</td><td>\code
-rot2 = rot1.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>
-
-\subsection TutorialGeoEulerAngles Euler angles
-<table class="tutorial_code">
-<tr><td style="max-width:30em;">
-Euler angles might be convenient to create rotation object.
-Since there exist 24 differents convensions, they are one
-the ahand 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>
-
-<a href="#" class="top">top</a>\section TutorialGeoTransformation Affine transformations
-In Eigen we have chosen to not distinghish between points and vectors such that all points are
-actually represented by displacement vectors from the origine ( \f$ \mathbf{p} \equiv \mathbf{p}-0 \f$ ).
-With that in mind, real points and vector distinguish when the rotation is applied.
-<table class="tutorial_code">
-<tr><td></td><td>\b 3D </td><td>\b 2D </td></tr>
-<tr><td>Creation \n <span class="note">rot2D can also be an angle in radian</span></td><td>\code
-Transform3f t;
-t.fromPositionOrientationScale(
- pos,any_3D_rotation,Vector3f(sx,sy,sz)); \endcode</td><td>\code
-Transform2f t;
-t.fromPositionOrientationScale(
- pos,any_2D_rotation,Vector2f(sx,sy)); \endcode</td></tr>
-<tr><td>Apply the transformation to a \b point </td><td>\code
-Vector3f p1, p2;
-p2 = t * p1;\endcode</td><td>\code
-Vector2f p1, p2;
-p2 = t * p1;\endcode</td></tr>
-<tr><td>Apply the transformation to a \b vector </td><td>\code
-Vector3f v1, v2;
-v2 = t.linear() * v1;\endcode</td><td>\code
-Vector2f v1, v2;
-v2 = t.linear() * v1;\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 colspan="2">\code
-Matrix{3,2}f normalMatrix = t.linear().inverse().transpose();
-n2 = (normalMatrix * n1).normalize();\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 colspan="2">\code
-n2 = t.linear() * n1;\endcode</td></tr>
-<tr><td>Concatenate two transformations</td><td colspan="2">\code
-t3 = t1 * t2;\endcode</td></tr>
-<tr><td>OpenGL compatibility</td><td>\code
-glLoadMatrixf(t.data());\endcode</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>
-<tr><td colspan="3">\b Component \b accessors</td></tr>
-<tr><td>full read-write access to the internal matrix</td><td>\code
-t.matrix() = mat4x4;
-mat4x4 = t.matrix();
-\endcode</td><td>\code
-t.matrix() = mat3x3;
-mat3x3 = t.matrix();
-\endcode</td></tr>
-<tr><td>coefficient accessors</td><td colspan="2">\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() = vec3;
-vec3 = t.translation();
-\endcode</td><td>\code
-t.translation() = vec2;
-vec2 = t.translation();
-\endcode</td></tr>
-<tr><td>linear part</td><td>\code
-t.linear() = mat3x3;
-mat3x3 = t.linear();
-\endcode</td><td>\code
-t.linear() = mat2x2;
-mat2x2 = t.linear();
-\endcode</td></tr>
-</table>
-
-\b Transformation \b creation \n
-Eigen's geometry module offer two different ways to build and update transformation objects.
-<table class="tutorial_code">
-<tr><td></td><td>\b procedurale \b API </td><td>\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, ...);
-t2 = t1 * Translation_(vec);
-t = Translation_(tx, ty, ...) * t;
-\endcode</td></tr>
-<tr><td>Rotation \n <span class="note">In 2D, any_rotation can also \n be an angle in radian</span></td><td>\code
-t.rotate(any_rotation);
-
-t.prerotate(any_rotation);
-\endcode</td><td>\code
-t *= any_rotation;
-t2 = t1 * 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, ...);
-t2 = t1 * Scaling_(vec);
-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>
-
-*/
-
-
-
-
-
/** \page TutorialAdvancedLinearAlgebra Tutorial 3/3 - Advanced linear algebra
\ingroup Tutorial
diff --git a/doc/TutorialGeometry.dox b/doc/TutorialGeometry.dox
index 69ca78997..8742341f5 100644
--- a/doc/TutorialGeometry.dox
+++ b/doc/TutorialGeometry.dox
@@ -23,13 +23,13 @@ namely 2D and 3D rotations and affine transformations.
<table class="tutorial_code">
<tr><td>Transformation type</td><td>Typical initialization code</td></tr>
<tr><td>
-2D rotation from an angle</td><td>\code
+\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 angle + axis</td><td>\code
+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 quaternion</td><td>\code
+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
diff --git a/doc/eigendoxy.css b/doc/eigendoxy.css
index 1e20503c5..cd7e97868 100644
--- a/doc/eigendoxy.css
+++ b/doc/eigendoxy.css
@@ -461,6 +461,23 @@ DIV.navigation {
padding-top : 5px;
}
+img {
+ border: 0;
+}
+table {
+ border-collapse: collapse;
+ border-style: none;
+ font-size: 1em;
+}
+th {
+ text-align: left;
+ padding-right: 1em;
+/* border: #cccccc dashed; */
+/* border-style: dashed; */
+/* border-width: 0 0 3px 0; */
+}
+
+
TABLE.noborder {
border-bottom-style : none;
border-left-style : none;
@@ -480,22 +497,18 @@ TABLE.noborder TD {
vertical-align: top;
}
-TABLE.tutorial_code {
- border-bottom-style : none;
- border-left-style : dashed;
- border-right-style : dashed;
- border-top-style : dashed ;
- border-spacing : 0px 0px;
+table.tutorial_code {
+ border-collapse: collapse;
empty-cells : show;
margin: 4pt 0 0 0;
padding: 0 0 0 0;
}
-TABLE.tutorial_code TD {
- border-bottom-style : dashed;
- border-left-style : none;
- border-right-style : none;
- border-top-style : none;
- border-spacing : 0px 0px;
+table.tutorial_code tr {
+ border-style: none dashed none dashed;
+ border-width: 1px;
+}
+table.tutorial_code td {
+ border-style: dashed none dashed none;
empty-cells : show;
margin: 0 0 0 0;
padding: 2pt 5pt 2pt 5pt;