mirror of
https://gitlab.com/libeigen/eigen.git
synced 2024-12-15 07:10:37 +08:00
* rename JacobiRotation => PlanarRotation
* move the makeJacobi and make_givens_* to PlanarRotation * rename applyJacobi* => apply*
This commit is contained in:
parent
496ea63972
commit
b83654b5d0
@ -116,6 +116,7 @@ inline float ei_imag(float) { return 0.f; }
|
||||
inline float ei_conj(float x) { return x; }
|
||||
inline float ei_abs(float x) { return std::abs(x); }
|
||||
inline float ei_abs2(float x) { return x*x; }
|
||||
inline float ei_norm1(float x) { return ei_abs(x); }
|
||||
inline float ei_sqrt(float x) { return std::sqrt(x); }
|
||||
inline float ei_exp(float x) { return std::exp(x); }
|
||||
inline float ei_log(float x) { return std::log(x); }
|
||||
@ -164,6 +165,7 @@ inline double ei_imag(double) { return 0.; }
|
||||
inline double ei_conj(double x) { return x; }
|
||||
inline double ei_abs(double x) { return std::abs(x); }
|
||||
inline double ei_abs2(double x) { return x*x; }
|
||||
inline double ei_norm1(double x) { return ei_abs(x); }
|
||||
inline double ei_sqrt(double x) { return std::sqrt(x); }
|
||||
inline double ei_exp(double x) { return std::exp(x); }
|
||||
inline double ei_log(double x) { return std::log(x); }
|
||||
@ -212,6 +214,7 @@ inline float& ei_imag_ref(std::complex<float>& x) { return reinterpret_cast<floa
|
||||
inline std::complex<float> ei_conj(const std::complex<float>& x) { return std::conj(x); }
|
||||
inline float ei_abs(const std::complex<float>& x) { return std::abs(x); }
|
||||
inline float ei_abs2(const std::complex<float>& x) { return std::norm(x); }
|
||||
inline float ei_norm1(const std::complex<float> &x) { return(ei_abs(x.real()) + ei_abs(x.imag())); }
|
||||
inline std::complex<float> ei_exp(std::complex<float> x) { return std::exp(x); }
|
||||
inline std::complex<float> ei_sin(std::complex<float> x) { return std::sin(x); }
|
||||
inline std::complex<float> ei_cos(std::complex<float> x) { return std::cos(x); }
|
||||
@ -248,6 +251,7 @@ inline double& ei_imag_ref(std::complex<double>& x) { return reinterpret_cast<do
|
||||
inline std::complex<double> ei_conj(const std::complex<double>& x) { return std::conj(x); }
|
||||
inline double ei_abs(const std::complex<double>& x) { return std::abs(x); }
|
||||
inline double ei_abs2(const std::complex<double>& x) { return std::norm(x); }
|
||||
inline double ei_norm1(const std::complex<double> &x) { return(ei_abs(x.real()) + ei_abs(x.imag())); }
|
||||
inline std::complex<double> ei_exp(std::complex<double> x) { return std::exp(x); }
|
||||
inline std::complex<double> ei_sin(std::complex<double> x) { return std::sin(x); }
|
||||
inline std::complex<double> ei_cos(std::complex<double> x) { return std::cos(x); }
|
||||
|
@ -803,11 +803,10 @@ template<typename Derived> class MatrixBase
|
||||
|
||||
///////// Jacobi module /////////
|
||||
|
||||
template<typename JacobiScalar>
|
||||
void applyJacobiOnTheLeft(int p, int q, const JacobiRotation<JacobiScalar>& j);
|
||||
template<typename JacobiScalar>
|
||||
void applyJacobiOnTheRight(int p, int q, const JacobiRotation<JacobiScalar>& j);
|
||||
bool makeJacobi(int p, int q, JacobiRotation<Scalar> *j) const;
|
||||
template<typename OtherScalar>
|
||||
void applyOnTheLeft(int p, int q, const PlanarRotation<OtherScalar>& j);
|
||||
template<typename OtherScalar>
|
||||
void applyOnTheRight(int p, int q, const PlanarRotation<OtherScalar>& j);
|
||||
|
||||
#ifdef EIGEN_MATRIXBASE_PLUGIN
|
||||
#include EIGEN_MATRIXBASE_PLUGIN
|
||||
|
@ -123,7 +123,7 @@ template<typename MatrixType> class SVD;
|
||||
template<typename MatrixType, unsigned int Options = 0> class JacobiSVD;
|
||||
template<typename MatrixType, int UpLo = LowerTriangular> class LLT;
|
||||
template<typename MatrixType> class LDLT;
|
||||
template<typename Scalar> class JacobiRotation;
|
||||
template<typename Scalar> class PlanarRotation;
|
||||
|
||||
// Geometry module:
|
||||
template<typename Derived, int _Dim> class RotationBase;
|
||||
|
@ -27,97 +27,72 @@
|
||||
#define EIGEN_JACOBI_H
|
||||
|
||||
/** \ingroup Jacobi
|
||||
* \class JacobiRotation
|
||||
* \class PlanarRotation
|
||||
* \brief Represents a rotation in the plane from a cosine-sine pair.
|
||||
*
|
||||
* This class represents a Jacobi rotation which is also known as a Givens rotation.
|
||||
* 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
|
||||
* 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::applyJacobiOnTheLeft(), MatrixBase::applyJacobiOnTheRight()
|
||||
* \sa MatrixBase::makeJacobi(), MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
|
||||
*/
|
||||
template<typename Scalar> class JacobiRotation
|
||||
template<typename Scalar> class PlanarRotation
|
||||
{
|
||||
public:
|
||||
/** Default constructor without any initialization. */
|
||||
JacobiRotation() {}
|
||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||
|
||||
/** Construct a Jacobi rotation from a cosine-sine pair (\a c, \c s). */
|
||||
JacobiRotation(const Scalar& c, const Scalar& s) : m_c(c), m_s(s) {}
|
||||
/** Default constructor without any initialization. */
|
||||
PlanarRotation() {}
|
||||
|
||||
/** Construct a planar rotation from a cosine-sine pair (\a c, \c s). */
|
||||
PlanarRotation(const Scalar& c, const Scalar& s) : m_c(c), m_s(s) {}
|
||||
|
||||
Scalar& c() { return m_c; }
|
||||
Scalar c() const { return m_c; }
|
||||
Scalar& s() { return m_s; }
|
||||
Scalar s() const { return m_s; }
|
||||
|
||||
/** Concatenates two Jacobi rotation */
|
||||
JacobiRotation operator*(const JacobiRotation& other)
|
||||
/** Concatenates two planar rotation */
|
||||
PlanarRotation operator*(const PlanarRotation& other)
|
||||
{
|
||||
return JacobiRotation(m_c * other.m_c - ei_conj(m_s) * other.m_s,
|
||||
return PlanarRotation(m_c * other.m_c - ei_conj(m_s) * other.m_s,
|
||||
ei_conj(m_c * ei_conj(other.m_s) + ei_conj(m_s) * ei_conj(other.m_c)));
|
||||
}
|
||||
|
||||
/** Returns the transposed transformation */
|
||||
JacobiRotation transpose() const { return JacobiRotation(m_c, -ei_conj(m_s)); }
|
||||
PlanarRotation transpose() const { return PlanarRotation(m_c, -ei_conj(m_s)); }
|
||||
|
||||
/** Returns the adjoint transformation */
|
||||
JacobiRotation adjoint() const { return JacobiRotation(ei_conj(m_c), -m_s); }
|
||||
PlanarRotation adjoint() const { return PlanarRotation(ei_conj(m_c), -m_s); }
|
||||
|
||||
template<typename Derived>
|
||||
bool makeJacobi(const MatrixBase<Derived>&, int p, int q);
|
||||
bool makeJacobi(RealScalar x, Scalar y, RealScalar z);
|
||||
|
||||
void makeGivens(const Scalar& p, const Scalar& q, Scalar* z=0);
|
||||
|
||||
protected:
|
||||
void makeGivens(const Scalar& p, const Scalar& q, Scalar* z, ei_meta_true);
|
||||
void makeGivens(const Scalar& p, const Scalar& q, Scalar* z, ei_meta_false);
|
||||
|
||||
Scalar m_c, m_s;
|
||||
};
|
||||
|
||||
/** 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::applyJacobiOnTheLeft(), MatrixBase::applyJacobiOnTheRight()
|
||||
*/
|
||||
template<typename VectorX, typename VectorY, typename JacobiScalar>
|
||||
void ei_apply_rotation_in_the_plane(VectorX& _x, VectorY& _y, const JacobiRotation<JacobiScalar>& 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,
|
||||
* with \f$ B = \left ( \begin{array}{cc} \text{*this.row}(p) \\ \text{*this.row}(q) \end{array} \right ) \f$.
|
||||
*
|
||||
* \sa class JacobiRotation, MatrixBase::applyJacobiOnTheRight(), ei_apply_rotation_in_the_plane()
|
||||
*/
|
||||
template<typename Derived>
|
||||
template<typename JacobiScalar>
|
||||
inline void MatrixBase<Derived>::applyJacobiOnTheLeft(int p, int q, const JacobiRotation<JacobiScalar>& j)
|
||||
{
|
||||
RowXpr x(row(p));
|
||||
RowXpr y(row(q));
|
||||
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
|
||||
* with \f$ B = \left ( \begin{array}{cc} \text{*this.col}(p) & \text{*this.col}(q) \end{array} \right ) \f$.
|
||||
*
|
||||
* \sa class JacobiRotation, MatrixBase::applyJacobiOnTheLeft(), ei_apply_rotation_in_the_plane()
|
||||
*/
|
||||
template<typename Derived>
|
||||
template<typename JacobiScalar>
|
||||
inline void MatrixBase<Derived>::applyJacobiOnTheRight(int p, int q, const JacobiRotation<JacobiScalar>& j)
|
||||
{
|
||||
ColXpr x(col(p));
|
||||
ColXpr y(col(q));
|
||||
ei_apply_rotation_in_the_plane(x, y, j.transpose());
|
||||
}
|
||||
|
||||
/** Computes the Jacobi rotation \a J such that applying \a J on both the right and left sides of the 2x2 matrix
|
||||
/** 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$
|
||||
*
|
||||
* \sa MatrixBase::makeJacobi(), MatrixBase::applyJacobiOnTheLeft(), MatrixBase::applyJacobiOnTheRight()
|
||||
* \sa MatrixBase::makeJacobi(), MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
|
||||
*/
|
||||
template<typename Scalar>
|
||||
bool ei_makeJacobi(typename NumTraits<Scalar>::Real x, Scalar y, typename NumTraits<Scalar>::Real z, JacobiRotation<Scalar> *j)
|
||||
bool PlanarRotation<Scalar>::makeJacobi(RealScalar x, Scalar y, RealScalar z)
|
||||
{
|
||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||
if(y == Scalar(0))
|
||||
{
|
||||
j->c() = Scalar(1);
|
||||
j->s() = Scalar(0);
|
||||
m_c = Scalar(1);
|
||||
m_s = Scalar(0);
|
||||
return false;
|
||||
}
|
||||
else
|
||||
@ -135,26 +110,132 @@ bool ei_makeJacobi(typename NumTraits<Scalar>::Real x, Scalar y, typename NumTra
|
||||
}
|
||||
RealScalar sign_t = t > 0 ? 1 : -1;
|
||||
RealScalar n = RealScalar(1) / ei_sqrt(ei_abs2(t)+1);
|
||||
j->s() = - sign_t * (ei_conj(y) / ei_abs(y)) * ei_abs(t) * n;
|
||||
j->c() = n;
|
||||
m_s = - sign_t * (ei_conj(y) / ei_abs(y)) * ei_abs(t) * n;
|
||||
m_c = n;
|
||||
return true;
|
||||
}
|
||||
}
|
||||
|
||||
/** Computes the Jacobi rotation \a 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 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$
|
||||
*
|
||||
* \sa MatrixBase::ei_make_jacobi(), MatrixBase::applyJacobiOnTheLeft(), MatrixBase::applyJacobiOnTheRight()
|
||||
* \sa PlanarRotation::makeJacobi(RealScalar, Scalar, RealScalar), MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
|
||||
*/
|
||||
template<typename Scalar>
|
||||
template<typename Derived>
|
||||
inline bool MatrixBase<Derived>::makeJacobi(int p, int q, JacobiRotation<Scalar> *j) const
|
||||
inline bool PlanarRotation<Scalar>::makeJacobi(const MatrixBase<Derived>& m, int p, int q)
|
||||
{
|
||||
return ei_makeJacobi(ei_real(coeff(p,p)), coeff(p,q), ei_real(coeff(q,q)), j);
|
||||
return makeJacobi(ei_real(m.coeff(p,p)), m.coeff(p,q), ei_real(m.coeff(q,q)));
|
||||
}
|
||||
|
||||
template<typename VectorX, typename VectorY, typename JacobiScalar>
|
||||
void /*EIGEN_DONT_INLINE*/ ei_apply_rotation_in_the_plane(VectorX& _x, VectorY& _y, const JacobiRotation<JacobiScalar>& j)
|
||||
/** Makes \c *this as a Givens rotation \c G such that applying \f$ G^* \f$ to the left of the vector
|
||||
* \f$ V = \left ( \begin{array}{c} p \\ q \end{array} \right )\f$ yields:
|
||||
* \f$ G^* V = \left ( \begin{array}{c} z \\ 0 \end{array} \right )\f$.
|
||||
*
|
||||
* 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.
|
||||
*
|
||||
* \sa MatrixBase::applyOnTheLeft(), MatrixBase::applyOnTheRight()
|
||||
*/
|
||||
template<typename Scalar>
|
||||
void PlanarRotation<Scalar>::makeGivens(const Scalar& p, const Scalar& q, Scalar* z)
|
||||
{
|
||||
makeGivens(p, q, z, typename ei_meta_if<NumTraits<Scalar>::IsComplex, ei_meta_true, ei_meta_false>::ret());
|
||||
}
|
||||
|
||||
|
||||
// specialization for complexes
|
||||
template<typename Scalar>
|
||||
void PlanarRotation<Scalar>::makeGivens(const Scalar& p, const Scalar& q, Scalar* z, ei_meta_true)
|
||||
{
|
||||
RealScalar scale, absx, absxy;
|
||||
if(q==Scalar(0))
|
||||
{
|
||||
// return identity
|
||||
m_c = Scalar(1);
|
||||
m_s = Scalar(0);
|
||||
if(z) *z = p;
|
||||
}
|
||||
else
|
||||
{
|
||||
scale = ei_norm1(p);
|
||||
absx = scale * ei_sqrt(ei_abs2(p/scale));
|
||||
scale = ei_abs(scale) + ei_norm1(q);
|
||||
absxy = scale * ei_sqrt((absx/scale)*(absx/scale) + ei_abs2(q/scale));
|
||||
m_c = Scalar(absx / absxy);
|
||||
Scalar np = p/absx;
|
||||
m_s = -ei_conj(np) * q / absxy;
|
||||
if(z) *z = np * absxy;
|
||||
}
|
||||
}
|
||||
|
||||
// specialization for reals
|
||||
template<typename Scalar>
|
||||
void PlanarRotation<Scalar>::makeGivens(const Scalar& p, const Scalar& q, Scalar* z, ei_meta_false)
|
||||
{
|
||||
// from Golub's "Matrix Computations", algorithm 5.1.3
|
||||
if(q==0)
|
||||
{
|
||||
m_c = 1; m_s = 0;
|
||||
}
|
||||
else if(ei_abs(q)>ei_abs(p))
|
||||
{
|
||||
Scalar t = -p/q;
|
||||
m_s = Scalar(1)/ei_sqrt(1+t*t);
|
||||
m_c = m_s * t;
|
||||
}
|
||||
else
|
||||
{
|
||||
Scalar t = -q/p;
|
||||
m_c = Scalar(1)/ei_sqrt(1+t*t);
|
||||
m_s = m_c * t;
|
||||
}
|
||||
}
|
||||
|
||||
/****************************************************************************************
|
||||
* 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:
|
||||
* \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()
|
||||
*/
|
||||
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,
|
||||
* 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()
|
||||
*/
|
||||
template<typename Derived>
|
||||
template<typename OtherScalar>
|
||||
inline void MatrixBase<Derived>::applyOnTheLeft(int p, int q, const PlanarRotation<OtherScalar>& j)
|
||||
{
|
||||
RowXpr x(row(p));
|
||||
RowXpr y(row(q));
|
||||
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
|
||||
* 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()
|
||||
*/
|
||||
template<typename Derived>
|
||||
template<typename OtherScalar>
|
||||
inline void MatrixBase<Derived>::applyOnTheRight(int p, int q, const PlanarRotation<OtherScalar>& j)
|
||||
{
|
||||
ColXpr x(col(p));
|
||||
ColXpr y(col(q));
|
||||
ei_apply_rotation_in_the_plane(x, y, j.transpose());
|
||||
}
|
||||
|
||||
|
||||
template<typename VectorX, typename VectorY, typename OtherScalar>
|
||||
void /*EIGEN_DONT_INLINE*/ ei_apply_rotation_in_the_plane(VectorX& _x, VectorY& _y, const PlanarRotation<OtherScalar>& j)
|
||||
{
|
||||
typedef typename VectorX::Scalar Scalar;
|
||||
ei_assert(_x.size() == _y.size());
|
||||
|
@ -80,34 +80,6 @@ template<typename _MatrixType> class ComplexSchur
|
||||
bool m_isInitialized;
|
||||
};
|
||||
|
||||
// computes the plane rotation G such that G' x |p| = | c s' |' |p| = |z|
|
||||
// |q| |-s c' | |q| |0|
|
||||
// and returns z if requested. Note that the returned c is real.
|
||||
template<typename T> void ei_make_givens(const std::complex<T>& p, const std::complex<T>& q,
|
||||
JacobiRotation<std::complex<T> >& rot, std::complex<T>* z=0)
|
||||
{
|
||||
typedef std::complex<T> Complex;
|
||||
T scale, absx, absxy;
|
||||
if(p==Complex(0))
|
||||
{
|
||||
// return identity
|
||||
rot.c() = Complex(1,0);
|
||||
rot.s() = Complex(0,0);
|
||||
if(z) *z = p;
|
||||
}
|
||||
else
|
||||
{
|
||||
scale = cnorm1(p);
|
||||
absx = scale * ei_sqrt(ei_abs2(p/scale));
|
||||
scale = ei_abs(scale) + cnorm1(q);
|
||||
absxy = scale * ei_sqrt((absx/scale)*(absx/scale) + ei_abs2(q/scale));
|
||||
rot.c() = Complex(absx / absxy);
|
||||
Complex np = p/absx;
|
||||
rot.s() = -ei_conj(np) * q / absxy;
|
||||
if(z) *z = np * absxy;
|
||||
}
|
||||
}
|
||||
|
||||
template<typename MatrixType>
|
||||
void ComplexSchur<MatrixType>::compute(const MatrixType& matrix)
|
||||
{
|
||||
@ -133,8 +105,8 @@ void ComplexSchur<MatrixType>::compute(const MatrixType& matrix)
|
||||
//locate the range in which to iterate
|
||||
while(iu > 0)
|
||||
{
|
||||
d = cnorm1(m_matT.coeffRef(iu,iu)) + cnorm1(m_matT.coeffRef(iu-1,iu-1));
|
||||
sd = cnorm1(m_matT.coeffRef(iu,iu-1));
|
||||
d = ei_norm1(m_matT.coeffRef(iu,iu)) + ei_norm1(m_matT.coeffRef(iu-1,iu-1));
|
||||
sd = ei_norm1(m_matT.coeffRef(iu,iu-1));
|
||||
|
||||
if(sd >= eps * d) break; // FIXME : precision criterion ??
|
||||
|
||||
@ -156,8 +128,8 @@ void ComplexSchur<MatrixType>::compute(const MatrixType& matrix)
|
||||
while( il > 0 )
|
||||
{
|
||||
// check if the current 2x2 block on the diagonal is upper triangular
|
||||
d = cnorm1(m_matT.coeffRef(il,il)) + cnorm1(m_matT.coeffRef(il-1,il-1));
|
||||
sd = cnorm1(m_matT.coeffRef(il,il-1));
|
||||
d = ei_norm1(m_matT.coeffRef(il,il)) + ei_norm1(m_matT.coeffRef(il-1,il-1));
|
||||
sd = ei_norm1(m_matT.coeffRef(il,il-1));
|
||||
|
||||
if(sd < eps * d) break; // FIXME : precision criterion ??
|
||||
|
||||
@ -179,32 +151,32 @@ void ComplexSchur<MatrixType>::compute(const MatrixType& matrix)
|
||||
r1 = (b+disc)/RealScalar(2);
|
||||
r2 = (b-disc)/RealScalar(2);
|
||||
|
||||
if(cnorm1(r1) > cnorm1(r2))
|
||||
if(ei_norm1(r1) > ei_norm1(r2))
|
||||
r2 = c/r1;
|
||||
else
|
||||
r1 = c/r2;
|
||||
|
||||
if(cnorm1(r1-t.coeff(1,1)) < cnorm1(r2-t.coeff(1,1)))
|
||||
if(ei_norm1(r1-t.coeff(1,1)) < ei_norm1(r2-t.coeff(1,1)))
|
||||
kappa = sf * r1;
|
||||
else
|
||||
kappa = sf * r2;
|
||||
|
||||
// perform the QR step using Givens rotations
|
||||
JacobiRotation<Complex> rot;
|
||||
ei_make_givens(m_matT.coeff(il,il) - kappa, m_matT.coeff(il+1,il), rot);
|
||||
PlanarRotation<Complex> rot;
|
||||
rot.makeGivens(m_matT.coeff(il,il) - kappa, m_matT.coeff(il+1,il));
|
||||
|
||||
for(int i=il ; i<iu ; i++)
|
||||
{
|
||||
m_matT.block(0,i,n,n-i).applyJacobiOnTheLeft(i, i+1, rot.adjoint());
|
||||
m_matT.block(0,0,std::min(i+2,iu)+1,n).applyJacobiOnTheRight(i, i+1, rot);
|
||||
m_matU.applyJacobiOnTheRight(i, i+1, rot);
|
||||
m_matT.block(0,i,n,n-i).applyOnTheLeft(i, i+1, rot.adjoint());
|
||||
m_matT.block(0,0,std::min(i+2,iu)+1,n).applyOnTheRight(i, i+1, rot);
|
||||
m_matU.applyOnTheRight(i, i+1, rot);
|
||||
|
||||
if(i != iu-1)
|
||||
{
|
||||
int i1 = i+1;
|
||||
int i2 = i+2;
|
||||
|
||||
ei_make_givens(m_matT.coeffRef(i1,i), m_matT.coeffRef(i2,i), rot, &m_matT.coeffRef(i1,i));
|
||||
rot.makeGivens(m_matT.coeffRef(i1,i), m_matT.coeffRef(i2,i), &m_matT.coeffRef(i1,i));
|
||||
m_matT.coeffRef(i2,i) = Complex(0);
|
||||
}
|
||||
}
|
||||
@ -223,13 +195,6 @@ void ComplexSchur<MatrixType>::compute(const MatrixType& matrix)
|
||||
m_isInitialized = true;
|
||||
}
|
||||
|
||||
// norm1 of complex numbers
|
||||
template<typename T>
|
||||
T cnorm1(const std::complex<T> &Z)
|
||||
{
|
||||
return(ei_abs(Z.real()) + ei_abs(Z.imag()));
|
||||
}
|
||||
|
||||
/**
|
||||
* Computes the principal value of the square root of the complex \a z.
|
||||
*/
|
||||
|
@ -135,28 +135,6 @@ template<typename _MatrixType> class SelfAdjointEigenSolver
|
||||
|
||||
#ifndef EIGEN_HIDE_HEAVY_CODE
|
||||
|
||||
// from Golub's "Matrix Computations", algorithm 5.1.3
|
||||
template<typename Scalar>
|
||||
static void ei_givens_rotation(Scalar a, Scalar b, Scalar& c, Scalar& s)
|
||||
{
|
||||
if (b==0)
|
||||
{
|
||||
c = 1; s = 0;
|
||||
}
|
||||
else if (ei_abs(b)>ei_abs(a))
|
||||
{
|
||||
Scalar t = -a/b;
|
||||
s = Scalar(1)/ei_sqrt(1+t*t);
|
||||
c = s * t;
|
||||
}
|
||||
else
|
||||
{
|
||||
Scalar t = -b/a;
|
||||
c = Scalar(1)/ei_sqrt(1+t*t);
|
||||
s = c * t;
|
||||
}
|
||||
}
|
||||
|
||||
/** \internal
|
||||
*
|
||||
* \qr_module
|
||||
@ -353,34 +331,33 @@ static void ei_tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, int st
|
||||
|
||||
for (int k = start; k < end; ++k)
|
||||
{
|
||||
RealScalar c, s;
|
||||
ei_givens_rotation(x, z, c, s);
|
||||
PlanarRotation<RealScalar> rot;
|
||||
rot.makeGivens(x, z);
|
||||
|
||||
// do T = G' T G
|
||||
RealScalar sdk = s * diag[k] + c * subdiag[k];
|
||||
RealScalar dkp1 = s * subdiag[k] + c * diag[k+1];
|
||||
RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
|
||||
RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k+1];
|
||||
|
||||
diag[k] = c * (c * diag[k] - s * subdiag[k]) - s * (c * subdiag[k] - s * diag[k+1]);
|
||||
diag[k+1] = s * sdk + c * dkp1;
|
||||
subdiag[k] = c * sdk - s * dkp1;
|
||||
diag[k] = rot.c() * (rot.c() * diag[k] - rot.s() * subdiag[k]) - rot.s() * (rot.c() * subdiag[k] - rot.s() * diag[k+1]);
|
||||
diag[k+1] = rot.s() * sdk + rot.c() * dkp1;
|
||||
subdiag[k] = rot.c() * sdk - rot.s() * dkp1;
|
||||
|
||||
if (k > start)
|
||||
subdiag[k - 1] = c * subdiag[k-1] - s * z;
|
||||
subdiag[k - 1] = rot.c() * subdiag[k-1] - rot.s() * z;
|
||||
|
||||
x = subdiag[k];
|
||||
|
||||
if (k < end - 1)
|
||||
{
|
||||
z = -s * subdiag[k+1];
|
||||
subdiag[k + 1] = c * subdiag[k+1];
|
||||
z = -rot.s() * subdiag[k+1];
|
||||
subdiag[k + 1] = rot.c() * subdiag[k+1];
|
||||
}
|
||||
|
||||
// apply the givens rotation to the unit matrix Q = Q * G
|
||||
// G only modifies the two columns k and k+1
|
||||
if (matrixQ)
|
||||
{
|
||||
Map<Matrix<Scalar,Dynamic,Dynamic> > q(matrixQ,n,n);
|
||||
q.applyJacobiOnTheRight(k,k+1,JacobiRotation<RealScalar>(c,s));
|
||||
q.applyOnTheRight(k,k+1,rot);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
@ -125,7 +125,7 @@ struct ei_svd_precondition_2x2_block_to_be_real<MatrixType, Options, true>
|
||||
static void run(MatrixType& work_matrix, JacobiSVD<MatrixType, Options>& svd, int p, int q)
|
||||
{
|
||||
Scalar z;
|
||||
JacobiRotation<Scalar> rot;
|
||||
PlanarRotation<Scalar> rot;
|
||||
RealScalar n = ei_sqrt(ei_abs2(work_matrix.coeff(p,p)) + ei_abs2(work_matrix.coeff(q,p)));
|
||||
if(n==0)
|
||||
{
|
||||
@ -140,8 +140,8 @@ struct ei_svd_precondition_2x2_block_to_be_real<MatrixType, Options, true>
|
||||
{
|
||||
rot.c() = ei_conj(work_matrix.coeff(p,p)) / n;
|
||||
rot.s() = work_matrix.coeff(q,p) / n;
|
||||
work_matrix.applyJacobiOnTheLeft(p,q,rot);
|
||||
if(ComputeU) svd.m_matrixU.applyJacobiOnTheRight(p,q,rot.adjoint());
|
||||
work_matrix.applyOnTheLeft(p,q,rot);
|
||||
if(ComputeU) svd.m_matrixU.applyOnTheRight(p,q,rot.adjoint());
|
||||
if(work_matrix.coeff(p,q) != Scalar(0))
|
||||
{
|
||||
Scalar z = ei_abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
|
||||
@ -160,13 +160,13 @@ struct ei_svd_precondition_2x2_block_to_be_real<MatrixType, Options, true>
|
||||
|
||||
template<typename MatrixType, typename RealScalar>
|
||||
void ei_real_2x2_jacobi_svd(const MatrixType& matrix, int p, int q,
|
||||
JacobiRotation<RealScalar> *j_left,
|
||||
JacobiRotation<RealScalar> *j_right)
|
||||
PlanarRotation<RealScalar> *j_left,
|
||||
PlanarRotation<RealScalar> *j_right)
|
||||
{
|
||||
Matrix<RealScalar,2,2> m;
|
||||
m << ei_real(matrix.coeff(p,p)), ei_real(matrix.coeff(p,q)),
|
||||
ei_real(matrix.coeff(q,p)), ei_real(matrix.coeff(q,q));
|
||||
JacobiRotation<RealScalar> rot1;
|
||||
PlanarRotation<RealScalar> rot1;
|
||||
RealScalar t = m.coeff(0,0) + m.coeff(1,1);
|
||||
RealScalar d = m.coeff(1,0) - m.coeff(0,1);
|
||||
if(t == RealScalar(0))
|
||||
@ -180,8 +180,8 @@ void ei_real_2x2_jacobi_svd(const MatrixType& matrix, int p, int q,
|
||||
rot1.c() = RealScalar(1) / ei_sqrt(1 + ei_abs2(u));
|
||||
rot1.s() = rot1.c() * u;
|
||||
}
|
||||
m.applyJacobiOnTheLeft(0,1,rot1);
|
||||
m.makeJacobi(0,1,j_right);
|
||||
m.applyOnTheLeft(0,1,rot1);
|
||||
j_right->makeJacobi(m,0,1);
|
||||
*j_left = rot1 * j_right->transpose();
|
||||
}
|
||||
|
||||
@ -232,14 +232,14 @@ sweep_again:
|
||||
{
|
||||
ei_svd_precondition_2x2_block_to_be_real<MatrixType, Options>::run(work_matrix, *this, p, q);
|
||||
|
||||
JacobiRotation<RealScalar> j_left, j_right;
|
||||
PlanarRotation<RealScalar> j_left, j_right;
|
||||
ei_real_2x2_jacobi_svd(work_matrix, p, q, &j_left, &j_right);
|
||||
|
||||
work_matrix.applyJacobiOnTheLeft(p,q,j_left);
|
||||
if(ComputeU) m_matrixU.applyJacobiOnTheRight(p,q,j_left.transpose());
|
||||
work_matrix.applyOnTheLeft(p,q,j_left);
|
||||
if(ComputeU) m_matrixU.applyOnTheRight(p,q,j_left.transpose());
|
||||
|
||||
work_matrix.applyJacobiOnTheRight(p,q,j_right);
|
||||
if(ComputeV) m_matrixV.applyJacobiOnTheRight(p,q,j_right);
|
||||
work_matrix.applyOnTheRight(p,q,j_right);
|
||||
if(ComputeV) m_matrixV.applyOnTheRight(p,q,j_right);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
@ -309,7 +309,7 @@ SVD<MatrixType>& SVD<MatrixType>::compute(const MatrixType& matrix)
|
||||
h = Scalar(1.0)/h;
|
||||
c = g*h;
|
||||
s = -f*h;
|
||||
V.applyJacobiOnTheRight(i,nm,JacobiRotation<Scalar>(c,s));
|
||||
V.applyOnTheRight(i,nm,PlanarRotation<Scalar>(c,s));
|
||||
}
|
||||
}
|
||||
z = W[k];
|
||||
@ -351,7 +351,7 @@ SVD<MatrixType>& SVD<MatrixType>::compute(const MatrixType& matrix)
|
||||
g = g*c - x*s;
|
||||
h = y*s;
|
||||
y *= c;
|
||||
V.applyJacobiOnTheRight(i,j,JacobiRotation<Scalar>(c,s));
|
||||
V.applyOnTheRight(i,j,PlanarRotation<Scalar>(c,s));
|
||||
|
||||
z = pythag(f,h);
|
||||
W[j] = z;
|
||||
@ -364,7 +364,7 @@ SVD<MatrixType>& SVD<MatrixType>::compute(const MatrixType& matrix)
|
||||
}
|
||||
f = c*g + s*y;
|
||||
x = c*y - s*g;
|
||||
A.applyJacobiOnTheRight(i,j,JacobiRotation<Scalar>(c,s));
|
||||
A.applyOnTheRight(i,j,PlanarRotation<Scalar>(c,s));
|
||||
}
|
||||
rv1[l] = 0.0;
|
||||
rv1[k] = f;
|
||||
|
Loading…
Reference in New Issue
Block a user