add examples for makeJacobi and makeGivens

This commit is contained in:
Gael Guennebaud 2009-09-03 11:17:16 +02:00
parent c893917d65
commit 16c7b1daab
6 changed files with 49 additions and 12 deletions

View File

@ -8,11 +8,15 @@
namespace Eigen {
/** \defgroup Jacobi_Module Jacobi module
* This module provides Jacobi rotations.
* This module provides Jacobi and Givens rotations.
*
* \code
* #include <Eigen/Jacobi>
* \endcode
*
* In addition to listed classes, it defines the two following MatrixBase methods to apply a Jacobi or Givens rotation:
* - MatrixBase::applyOnTheLeft()
* - MatrixBase::applyOnTheRight().
*/
#include "src/Jacobi/Jacobi.h"

View File

@ -26,16 +26,23 @@
#ifndef EIGEN_JACOBI_H
#define EIGEN_JACOBI_H
/** \ingroup Jacobi
/** \ingroup Jacobi_Module
* \jacobi_module
* \class PlanarRotation
* \brief Represents a rotation in the plane from a cosine-sine pair.
*
* This class represents a Jacobi or Givens rotation.
* This is a 2D clock-wise rotation in the plane \c J of angle \f$ \theta \f$ defined by
* This is a 2D rotation in the plane \c J of angle \f$ \theta \f$ defined by
* its cosine \c c and sine \c s as follow:
* \f$ J = \left ( \begin{array}{cc} c & \overline s \\ -s & \overline c \end{array} \right ) \f$
*
* \sa MatrixBase::makeJacobi(), MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
* You can apply the respective counter-clockwise rotation to a column vector \c v by
* applying its adjoint on the left: \f$ v = J^* v \f$ that translates to the following Eigen code:
* \code
* v.applyOnTheLeft(J.adjoint());
* \endcode
*
* \sa MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
*/
template<typename Scalar> class PlanarRotation
{
@ -79,11 +86,10 @@ template<typename Scalar> class PlanarRotation
Scalar m_c, m_s;
};
/** Makes \c *this as a Jacobi rotation \a J such that applying \a J on both the right and left sides of the 2x2 matrix
* \f$ B = \left ( \begin{array}{cc} x & y \\ * & z \end{array} \right )\f$ yields
* a diagonal matrix \f$ A = J^* B J \f$
/** Makes \c *this as a Jacobi rotation \a J such that applying \a J on both the right and left sides of the selfadjoint 2x2 matrix
* \f$ B = \left ( \begin{array}{cc} x & y \\ * & z \end{array} \right )\f$ yields a diagonal matrix \f$ A = J^* B J \f$
*
* \sa MatrixBase::makeJacobi(), MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
* \sa MatrixBase::makeJacobi(const MatrixBase<Derived>&, int, int), MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
*/
template<typename Scalar>
bool PlanarRotation<Scalar>::makeJacobi(RealScalar x, Scalar y, RealScalar z)
@ -116,10 +122,13 @@ bool PlanarRotation<Scalar>::makeJacobi(RealScalar x, Scalar y, RealScalar z)
}
}
/** Makes \c *this as a Jacobi rotation \c J such that applying \a J on both the right and left sides of the 2x2 matrix
/** Makes \c *this as a Jacobi rotation \c J such that applying \a J on both the right and left sides of the 2x2 selfadjoint matrix
* \f$ B = \left ( \begin{array}{cc} \text{this}_{pp} & \text{this}_{pq} \\ * & \text{this}_{qq} \end{array} \right )\f$ yields
* a diagonal matrix \f$ A = J^* B J \f$
*
* Example: \include Jacobi_makeJacobi.cpp
* Output: \verbinclude Jacobi_makeJacobi.out
*
* \sa PlanarRotation::makeJacobi(RealScalar, Scalar, RealScalar), MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
*/
template<typename Scalar>
@ -136,6 +145,9 @@ inline bool PlanarRotation<Scalar>::makeJacobi(const MatrixBase<Derived>& m, int
* The value of \a z is returned if \a z is not null (the default is null).
* Also note that G is built such that the cosine is always real.
*
* Example: \include Jacobi_makeGivens.cpp
* Output: \verbinclude Jacobi_makeGivens.out
*
* \sa MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
*/
template<typename Scalar>
@ -171,9 +183,11 @@ void PlanarRotation<Scalar>::makeGivens(const Scalar& p, const Scalar& q, Scalar
}
// specialization for reals
// TODO compute z
template<typename Scalar>
void PlanarRotation<Scalar>::makeGivens(const Scalar& p, const Scalar& q, Scalar* z, ei_meta_false)
{
ei_assert(z==0 && "not implemented yet");
// from Golub's "Matrix Computations", algorithm 5.1.3
if(q==0)
{
@ -197,7 +211,8 @@ void PlanarRotation<Scalar>::makeGivens(const Scalar& p, const Scalar& q, Scalar
* Implementation of MatrixBase methods
****************************************************************************************/
/** Applies the clock wise 2D rotation \a j to the set of 2D vectors of cordinates \a x and \a y:
/** \jacobi_module
* Applies the clock wise 2D rotation \a j to the set of 2D vectors of cordinates \a x and \a y:
* \f$ \left ( \begin{array}{cc} x \\ y \end{array} \right ) = J \left ( \begin{array}{cc} x \\ y \end{array} \right ) \f$
*
* \sa MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
@ -205,7 +220,8 @@ void PlanarRotation<Scalar>::makeGivens(const Scalar& p, const Scalar& q, Scalar
template<typename VectorX, typename VectorY, typename OtherScalar>
void ei_apply_rotation_in_the_plane(VectorX& _x, VectorY& _y, const PlanarRotation<OtherScalar>& j);
/** Applies the rotation in the plane \a j to the rows \a p and \a q of \c *this, i.e., it computes B = J * B,
/** \jacobi_module
* Applies the rotation in the plane \a j to the rows \a p and \a q of \c *this, i.e., it computes B = J * B,
* with \f$ B = \left ( \begin{array}{cc} \text{*this.row}(p) \\ \text{*this.row}(q) \end{array} \right ) \f$.
*
* \sa class PlanarRotation, MatrixBase::applyOnTheRight(), ei_apply_rotation_in_the_plane()
@ -219,7 +235,8 @@ inline void MatrixBase<Derived>::applyOnTheLeft(int p, int q, const PlanarRotati
ei_apply_rotation_in_the_plane(x, y, j);
}
/** Applies the rotation in the plane \a j to the columns \a p and \a q of \c *this, i.e., it computes B = B * J
/** \ingroup Jacobi_Module
* Applies the rotation in the plane \a j to the columns \a p and \a q of \c *this, i.e., it computes B = B * J
* with \f$ B = \left ( \begin{array}{cc} \text{*this.col}(p) & \text{*this.col}(q) \end{array} \right ) \f$.
*
* \sa class PlanarRotation, MatrixBase::applyOnTheLeft(), ei_apply_rotation_in_the_plane()

View File

@ -208,6 +208,7 @@ ALIASES = "only_for_vectors=This is only for vectors (either row-
"lu_module=This is defined in the %LU module. \code #include <Eigen/LU> \endcode" \
"cholesky_module=This is defined in the %Cholesky module. \code #include <Eigen/Cholesky> \endcode" \
"qr_module=This is defined in the %QR module. \code #include <Eigen/QR> \endcode" \
"jacobi_module=This is defined in the %Jacobi module. \code #include <Eigen/Jacobi> \endcode" \
"svd_module=This is defined in the %SVD module. \code #include <Eigen/SVD> \endcode" \
"geometry_module=This is defined in the %Geometry module. \code #include <Eigen/Geometry> \endcode" \
"leastsquares_module=This is defined in the %LeastSquares module. \code #include <Eigen/LeastSquares> \endcode" \

View File

@ -0,0 +1,6 @@
Vector2f v = Vector2f::Random();
PlanarRotation<float> G;
G.makeGivens(v.x(), v.y());
cout << "Here is the vector v:" << endl << v << endl;
v.applyOnTheLeft(0, 1, G.adjoint());
cout << "Here is the vector J' * v:" << endl << v << endl;

View File

@ -0,0 +1,8 @@
Matrix2f m = Matrix2f::Random();
m = (m + m.adjoint()).eval();
PlanarRotation<float> J;
J.makeJacobi(m, 0, 1);
cout << "Here is the matrix m:" << endl << m << endl;
m.applyOnTheLeft(0, 1, J.adjoint());
m.applyOnTheRight(0, 1, J);
cout << "Here is the matrix J' * m * J:" << endl << m << endl;

View File

@ -4,6 +4,7 @@
#include <Eigen/QR>
#include <Eigen/Cholesky>
#include <Eigen/Geometry>
#include <Eigen/Jacobi>
using namespace Eigen;
using namespace std;