MOAB
4.9.3pre
|
Rotation given by a cosine-sine pair. More...
#include <Jacobi.h>
Public Types | |
typedef NumTraits< Scalar >::Real | RealScalar |
Public Member Functions | |
JacobiRotation () | |
JacobiRotation (const Scalar &c, const Scalar &s) | |
Scalar & | c () |
Scalar | c () const |
Scalar & | s () |
Scalar | s () const |
JacobiRotation | operator* (const JacobiRotation &other) |
JacobiRotation | transpose () const |
JacobiRotation | adjoint () const |
template<typename Derived > | |
bool | makeJacobi (const MatrixBase< Derived > &, Index p, Index q) |
bool | makeJacobi (const RealScalar &x, const Scalar &y, const RealScalar &z) |
void | makeGivens (const Scalar &p, const Scalar &q, Scalar *z=0) |
Protected Member Functions | |
void | makeGivens (const Scalar &p, const Scalar &q, Scalar *z, internal::true_type) |
void | makeGivens (const Scalar &p, const Scalar &q, Scalar *z, internal::false_type) |
Protected Attributes | |
Scalar | m_c |
Scalar | m_s |
Rotation given by a cosine-sine pair.
This class represents a Jacobi or Givens rotation. This is a 2D rotation in the plane J
of angle defined by its cosine
c
and sine s
as follow:
You can apply the respective counter-clockwise rotation to a column vector v
by applying its adjoint on the left: that translates to the following Eigen code:
v.applyOnTheLeft(J.adjoint());
typedef NumTraits<Scalar>::Real Eigen::JacobiRotation< Scalar >::RealScalar |
Eigen::JacobiRotation< Scalar >::JacobiRotation | ( | ) | [inline] |
Eigen::JacobiRotation< Scalar >::JacobiRotation | ( | const Scalar & | c, |
const Scalar & | s | ||
) | [inline] |
JacobiRotation Eigen::JacobiRotation< Scalar >::adjoint | ( | ) | const [inline] |
Returns the adjoint transformation
Definition at line 62 of file Jacobi.h.
{ using numext::conj; return JacobiRotation(conj(m_c), -m_s); }
Scalar& Eigen::JacobiRotation< Scalar >::c | ( | ) | [inline] |
Scalar Eigen::JacobiRotation< Scalar >::c | ( | ) | const [inline] |
void Eigen::JacobiRotation< Scalar >::makeGivens | ( | const Scalar & | p, |
const Scalar & | q, | ||
Scalar * | z = 0 |
||
) |
Makes *this
as a Givens rotation G
such that applying to the left of the vector
yields:
.
The value of z is returned if z is not null (the default is null). Also note that G is built such that the cosine is always real.
Example:
Output:
This function implements the continuous Givens rotation generation algorithm found in Anderson (2000), Discontinuous Plane Rotations and the Symmetric Eigenvalue Problem. LAPACK Working Note 150, University of Tennessee, UT-CS-00-454, December 4, 2000.
Definition at line 148 of file Jacobi.h.
{ makeGivens(p, q, z, typename internal::conditional<NumTraits<Scalar>::IsComplex, internal::true_type, internal::false_type>::type()); }
void Eigen::JacobiRotation< Scalar >::makeGivens | ( | const Scalar & | p, |
const Scalar & | q, | ||
Scalar * | z, | ||
internal::true_type | |||
) | [protected] |
Definition at line 156 of file Jacobi.h.
{ using std::sqrt; using std::abs; using numext::conj; if(q==Scalar(0)) { m_c = numext::real(p)<0 ? Scalar(-1) : Scalar(1); m_s = 0; if(r) *r = m_c * p; } else if(p==Scalar(0)) { m_c = 0; m_s = -q/abs(q); if(r) *r = abs(q); } else { RealScalar p1 = numext::norm1(p); RealScalar q1 = numext::norm1(q); if(p1>=q1) { Scalar ps = p / p1; RealScalar p2 = numext::abs2(ps); Scalar qs = q / p1; RealScalar q2 = numext::abs2(qs); RealScalar u = sqrt(RealScalar(1) + q2/p2); if(numext::real(p)<RealScalar(0)) u = -u; m_c = Scalar(1)/u; m_s = -qs*conj(ps)*(m_c/p2); if(r) *r = p * u; } else { Scalar ps = p / q1; RealScalar p2 = numext::abs2(ps); Scalar qs = q / q1; RealScalar q2 = numext::abs2(qs); RealScalar u = q1 * sqrt(p2 + q2); if(numext::real(p)<RealScalar(0)) u = -u; p1 = abs(p); ps = p/p1; m_c = p1/u; m_s = -conj(ps) * (q/u); if(r) *r = ps * u; } } }
void Eigen::JacobiRotation< Scalar >::makeGivens | ( | const Scalar & | p, |
const Scalar & | q, | ||
Scalar * | z, | ||
internal::false_type | |||
) | [protected] |
Definition at line 215 of file Jacobi.h.
{ using std::sqrt; using std::abs; if(q==Scalar(0)) { m_c = p<Scalar(0) ? Scalar(-1) : Scalar(1); m_s = Scalar(0); if(r) *r = abs(p); } else if(p==Scalar(0)) { m_c = Scalar(0); m_s = q<Scalar(0) ? Scalar(1) : Scalar(-1); if(r) *r = abs(q); } else if(abs(p) > abs(q)) { Scalar t = q/p; Scalar u = sqrt(Scalar(1) + numext::abs2(t)); if(p<Scalar(0)) u = -u; m_c = Scalar(1)/u; m_s = -t * m_c; if(r) *r = p * u; } else { Scalar t = p/q; Scalar u = sqrt(Scalar(1) + numext::abs2(t)); if(q<Scalar(0)) u = -u; m_s = -Scalar(1)/u; m_c = -t * m_s; if(r) *r = q * u; } }
bool Eigen::JacobiRotation< Scalar >::makeJacobi | ( | const MatrixBase< Derived > & | m, |
Index | p, | ||
Index | q | ||
) | [inline] |
Makes *this
as a Jacobi rotation J
such that applying J on both the right and left sides of the 2x2 selfadjoint matrix yields a diagonal matrix
Example:
Output:
Definition at line 126 of file Jacobi.h.
{ return makeJacobi(numext::real(m.coeff(p,p)), m.coeff(p,q), numext::real(m.coeff(q,q))); }
bool Eigen::JacobiRotation< Scalar >::makeJacobi | ( | const RealScalar & | x, |
const Scalar & | y, | ||
const RealScalar & | z | ||
) |
Makes *this
as a Jacobi rotation J such that applying J on both the right and left sides of the selfadjoint 2x2 matrix yields a diagonal matrix
Definition at line 83 of file Jacobi.h.
{ using std::sqrt; using std::abs; typedef typename NumTraits<Scalar>::Real RealScalar; if(y == Scalar(0)) { m_c = Scalar(1); m_s = Scalar(0); return false; } else { RealScalar tau = (x-z)/(RealScalar(2)*abs(y)); RealScalar w = sqrt(numext::abs2(tau) + RealScalar(1)); RealScalar t; if(tau>RealScalar(0)) { t = RealScalar(1) / (tau + w); } else { t = RealScalar(1) / (tau - w); } RealScalar sign_t = t > RealScalar(0) ? RealScalar(1) : RealScalar(-1); RealScalar n = RealScalar(1) / sqrt(numext::abs2(t)+RealScalar(1)); m_s = - sign_t * (numext::conj(y) / abs(y)) * abs(t) * n; m_c = n; return true; } }
JacobiRotation Eigen::JacobiRotation< Scalar >::operator* | ( | const JacobiRotation< Scalar > & | other | ) | [inline] |
Scalar& Eigen::JacobiRotation< Scalar >::s | ( | ) | [inline] |
Scalar Eigen::JacobiRotation< Scalar >::s | ( | ) | const [inline] |
JacobiRotation Eigen::JacobiRotation< Scalar >::transpose | ( | ) | const [inline] |
Returns the transposed transformation
Definition at line 59 of file Jacobi.h.
{ using numext::conj; return JacobiRotation(m_c, -conj(m_s)); }
Scalar Eigen::JacobiRotation< Scalar >::m_c [protected] |
Scalar Eigen::JacobiRotation< Scalar >::m_s [protected] |