RealSchur: split computation in smaller functions.

This commit is contained in:
Jitse Niesen 2010-04-06 17:45:46 +01:00
parent 7dc56b3dfb
commit 9fad1e392b

View File

@ -93,7 +93,12 @@ template<typename _MatrixType> class RealSchur
EigenvalueType m_eivalues; EigenvalueType m_eivalues;
bool m_isInitialized; bool m_isInitialized;
void hqr2(); Scalar computeNormOfT();
int findSmallSubdiagEntry(int n, Scalar norm);
void computeShift(Scalar& x, Scalar& y, Scalar& w, int l, int n, Scalar& exshift, int iter);
void findTwoSmallSubdiagEntries(Scalar x, Scalar y, Scalar w, int l, int& m, int n, Scalar& p, Scalar& q, Scalar& r);
void doFrancisStep(int l, int m, int n, Scalar p, Scalar q, Scalar r, Scalar x, Scalar* workspace);
void splitOffTwoRows(int n, Scalar exshift);
}; };
@ -102,52 +107,71 @@ void RealSchur<MatrixType>::compute(const MatrixType& matrix)
{ {
assert(matrix.cols() == matrix.rows()); assert(matrix.cols() == matrix.rows());
// Reduce to Hessenberg form // Step 1. Reduce to Hessenberg form
// TODO skip Q if skipU = true // TODO skip Q if skipU = true
HessenbergDecomposition<MatrixType> hess(matrix); HessenbergDecomposition<MatrixType> hess(matrix);
m_matT = hess.matrixH(); m_matT = hess.matrixH();
m_matU = hess.matrixQ(); m_matU = hess.matrixQ();
// Reduce to Real Schur form // Step 2. Reduce to real Schur form
hqr2(); typedef Matrix<Scalar, ColsAtCompileTime, 1, Options, MaxColsAtCompileTime, 1> ColumnVectorType;
ColumnVectorType workspaceVector(m_matU.cols());
Scalar* workspace = &workspaceVector.coeffRef(0);
int n = m_matU.cols() - 1;
Scalar exshift = 0.0;
Scalar norm = computeNormOfT();
int iter = 0;
while (n >= 0)
{
int l = findSmallSubdiagEntry(n, norm);
// Check for convergence
if (l == n) // One root found
{
m_matT.coeffRef(n,n) = m_matT.coeff(n,n) + exshift;
m_eivalues.coeffRef(n) = ComplexScalar(m_matT.coeff(n,n), 0.0);
n--;
iter = 0;
}
else if (l == n-1) // Two roots found
{
splitOffTwoRows(n, exshift);
n = n - 2;
iter = 0;
}
else // No convergence yet
{
Scalar p = 0, q = 0, r = 0, x, y, w;
computeShift(x, y, w, l, n, exshift, iter);
iter = iter + 1; // (Could check iteration count here.)
int m;
findTwoSmallSubdiagEntries(x, y, w, l, m, n, p, q, r);
doFrancisStep(l, m, n, p, q, r, x, workspace);
} // check convergence
} // while (n >= 0)
m_isInitialized = true; m_isInitialized = true;
} }
// Compute matrix norm
template<typename MatrixType> template<typename MatrixType>
void RealSchur<MatrixType>::hqr2() inline typename MatrixType::Scalar RealSchur<MatrixType>::computeNormOfT()
{ {
typedef Matrix<Scalar, ColsAtCompileTime, 1, Options, MaxColsAtCompileTime, 1> ColumnVectorType;
// This is derived from the Algol procedure hqr2,
// by Martin and Wilkinson, Handbook for Auto. Comp.,
// Vol.ii-Linear Algebra, and the corresponding
// Fortran subroutine in EISPACK.
// Initialize
const int size = m_matU.cols(); const int size = m_matU.cols();
int n = size-1;
Scalar exshift = 0.0;
Scalar p=0, q=0, r=0;
ColumnVectorType workspaceVector(size);
Scalar* workspace = &workspaceVector.coeffRef(0);
// Compute matrix norm
// FIXME to be efficient the following would requires a triangular reduxion code // FIXME to be efficient the following would requires a triangular reduxion code
// Scalar norm = m_matT.upper().cwiseAbs().sum() + m_matT.corner(BottomLeft,n,n).diagonal().cwiseAbs().sum(); // Scalar norm = m_matT.upper().cwiseAbs().sum() + m_matT.corner(BottomLeft,size-1,size-1).diagonal().cwiseAbs().sum();
Scalar norm = 0.0; Scalar norm = 0.0;
for (int j = 0; j < size; ++j) for (int j = 0; j < size; ++j)
{
norm += m_matT.row(j).segment(std::max(j-1,0), size-std::max(j-1,0)).cwiseAbs().sum(); norm += m_matT.row(j).segment(std::max(j-1,0), size-std::max(j-1,0)).cwiseAbs().sum();
} return norm;
}
// Outer loop over eigenvalue index // Look for single small sub-diagonal element
int iter = 0; template<typename MatrixType>
while (n >= 0) inline int RealSchur<MatrixType>::findSmallSubdiagEntry(int n, Scalar norm)
{ {
// Look for single small sub-diagonal element
int l = n; int l = n;
while (l > 0) while (l > 0)
{ {
@ -158,21 +182,16 @@ void RealSchur<MatrixType>::hqr2()
break; break;
l--; l--;
} }
return l;
}
// Check for convergence template<typename MatrixType>
// One root found inline void RealSchur<MatrixType>::splitOffTwoRows(int n, Scalar exshift)
if (l == n) {
{ const int size = m_matU.cols();
m_matT.coeffRef(n,n) = m_matT.coeff(n,n) + exshift;
m_eivalues.coeffRef(n) = ComplexScalar(m_matT.coeff(n,n), 0.0);
n--;
iter = 0;
}
else if (l == n-1) // Two roots found
{
Scalar w = m_matT.coeff(n,n-1) * m_matT.coeff(n-1,n); Scalar w = m_matT.coeff(n,n-1) * m_matT.coeff(n-1,n);
p = (m_matT.coeff(n-1,n-1) - m_matT.coeff(n,n)) * Scalar(0.5); Scalar p = (m_matT.coeff(n-1,n-1) - m_matT.coeff(n,n)) * Scalar(0.5);
q = p * p + w; Scalar q = p * p + w;
Scalar z = ei_sqrt(ei_abs(q)); Scalar z = ei_sqrt(ei_abs(q));
m_matT.coeffRef(n,n) = m_matT.coeff(n,n) + exshift; m_matT.coeffRef(n,n) = m_matT.coeff(n,n) + exshift;
m_matT.coeffRef(n-1,n-1) = m_matT.coeff(n-1,n-1) + exshift; m_matT.coeffRef(n-1,n-1) = m_matT.coeff(n-1,n-1) + exshift;
@ -200,15 +219,15 @@ void RealSchur<MatrixType>::hqr2()
m_eivalues.coeffRef(n-1) = ComplexScalar(x + p, z); m_eivalues.coeffRef(n-1) = ComplexScalar(x + p, z);
m_eivalues.coeffRef(n) = ComplexScalar(x + p, -z); m_eivalues.coeffRef(n) = ComplexScalar(x + p, -z);
} }
n = n - 2; }
iter = 0;
} // Form shift
else // No convergence yet template<typename MatrixType>
{ inline void RealSchur<MatrixType>::computeShift(Scalar& x, Scalar& y, Scalar& w, int l, int n, Scalar& exshift, int iter)
// Form shift {
Scalar x = m_matT.coeff(n,n); x = m_matT.coeff(n,n);
Scalar y = 0.0; y = 0.0;
Scalar w = 0.0; w = 0.0;
if (l < n) if (l < n)
{ {
y = m_matT.coeff(n-1,n-1); y = m_matT.coeff(n-1,n-1);
@ -243,11 +262,13 @@ void RealSchur<MatrixType>::hqr2()
x = y = w = Scalar(0.964); x = y = w = Scalar(0.964);
} }
} }
}
iter = iter + 1; // (Could check iteration count here.) // Look for two consecutive small sub-diagonal elements
template<typename MatrixType>
// Look for two consecutive small sub-diagonal elements inline void RealSchur<MatrixType>::findTwoSmallSubdiagEntries(Scalar x, Scalar y, Scalar w, int l, int& m, int n, Scalar& p, Scalar& q, Scalar& r)
int m = n-2; {
m = n-2;
while (m >= l) while (m >= l)
{ {
Scalar z = m_matT.coeff(m,m); Scalar z = m_matT.coeff(m,m);
@ -278,8 +299,14 @@ void RealSchur<MatrixType>::hqr2()
if (i > m+2) if (i > m+2)
m_matT.coeffRef(i,i-3) = 0.0; m_matT.coeffRef(i,i-3) = 0.0;
} }
}
// Double QR step involving rows l:n and columns m:n
template<typename MatrixType>
inline void RealSchur<MatrixType>::doFrancisStep(int l, int m, int n, Scalar p, Scalar q, Scalar r, Scalar x, Scalar* workspace)
{
const int size = m_matU.cols();
// Double QR step involving rows l:n and columns m:n
for (int k = m; k <= n-1; ++k) for (int k = m; k <= n-1; ++k)
{ {
int notlast = (k != n-1); int notlast = (k != n-1);
@ -328,11 +355,8 @@ void RealSchur<MatrixType>::hqr2()
m_matT.block(0, k, std::min(n,k+3) + 1, 2).applyHouseholderOnTheRight(ess, p/s, workspace); m_matT.block(0, k, std::min(n,k+3) + 1, 2).applyHouseholderOnTheRight(ess, p/s, workspace);
m_matU.block(0, k, size, 2).applyHouseholderOnTheRight(ess, p/s, workspace); m_matU.block(0, k, size, 2).applyHouseholderOnTheRight(ess, p/s, workspace);
} }
} // (s != 0) } // (s != 0)
} // k loop } // k loop
} // check convergence
} // while (n >= 0)
} }
#endif // EIGEN_REAL_SCHUR_H #endif // EIGEN_REAL_SCHUR_H