// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2010 Benoit Jacob // // This Source Code Form is subject to the terms of the Mozilla // Public License v. 2.0. If a copy of the MPL was not distributed // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. #ifndef EIGEN_BIDIAGONALIZATION_H #define EIGEN_BIDIAGONALIZATION_H namespace Eigen { namespace internal { // UpperBidiagonalization will probably be replaced by a Bidiagonalization class, don't want to make it stable API. // At the same time, it's useful to keep for now as it's about the only thing that is testing the BandMatrix class. template class UpperBidiagonalization { public: typedef _MatrixType MatrixType; enum { RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, ColsAtCompileTimeMinusOne = internal::decrement_size::ret }; typedef typename MatrixType::Scalar Scalar; typedef typename MatrixType::RealScalar RealScalar; typedef typename MatrixType::Index Index; typedef Matrix RowVectorType; typedef Matrix ColVectorType; typedef BandMatrix BidiagonalType; typedef Matrix DiagVectorType; typedef Matrix SuperDiagVectorType; typedef HouseholderSequence< const MatrixType, CwiseUnaryOp, const Diagonal > > HouseholderUSequenceType; typedef HouseholderSequence< const typename internal::remove_all::type, Diagonal, OnTheRight > HouseholderVSequenceType; /** * \brief Default Constructor. * * The default constructor is useful in cases in which the user intends to * perform decompositions via Bidiagonalization::compute(const MatrixType&). */ UpperBidiagonalization() : m_householder(), m_bidiagonal(), m_isInitialized(false) {} UpperBidiagonalization(const MatrixType& matrix) : m_householder(matrix.rows(), matrix.cols()), m_bidiagonal(matrix.cols(), matrix.cols()), m_isInitialized(false) { compute(matrix); } UpperBidiagonalization& compute(const MatrixType& matrix); UpperBidiagonalization& computeUnblocked(const MatrixType& matrix); const MatrixType& householder() const { return m_householder; } const BidiagonalType& bidiagonal() const { return m_bidiagonal; } const HouseholderUSequenceType householderU() const { eigen_assert(m_isInitialized && "UpperBidiagonalization is not initialized."); return HouseholderUSequenceType(m_householder, m_householder.diagonal().conjugate()); } const HouseholderVSequenceType householderV() // const here gives nasty errors and i'm lazy { eigen_assert(m_isInitialized && "UpperBidiagonalization is not initialized."); return HouseholderVSequenceType(m_householder.conjugate(), m_householder.const_derived().template diagonal<1>()) .setLength(m_householder.cols()-1) .setShift(1); } protected: MatrixType m_householder; BidiagonalType m_bidiagonal; bool m_isInitialized; }; /** \internal * Helper routine for the block bidiagonal reduction. * This function reduces to bidiagonal form the left \c rows x \a blockSize vertical * and \a blockSize x \c cols horizontal panels of the \a A. The bottom-right block * is left unchanged, and the Arices X and Y. */ template void upperbidiagonalization_inplace_helper(MatrixType& A, typename MatrixType::RealScalar *diagonal, typename MatrixType::RealScalar *upper_diagonal, typename MatrixType::Index bs, Ref > X, Ref > Y) { typedef typename MatrixType::Index Index; typedef typename MatrixType::Scalar Scalar; typedef Ref > SubColumnType; typedef Ref, 0, InnerStride<> > SubRowType; typedef Ref > SubMatType; Index brows = A.rows(); Index bcols = A.cols(); Scalar tau_u, tau_u_prev(0), tau_v; for(Index k = 0; k < bs; ++k) { Index remainingRows = brows - k; Index remainingSizeInBlock = bs - k - 1; Index remainingCols = bcols - k - 1; SubMatType X_k1( X.block(k,0, remainingRows,k) ); SubMatType V_k1( A.block(k,0, remainingRows,k) ); // 1 - update the k-th column of A SubColumnType v_k = A.col(k).tail(remainingRows); v_k -= V_k1 * Y.row(k).head(k).adjoint(); if(k) v_k -= X_k1 * A.col(k).head(k); // 2 - construct left Householder transform in-place v_k.makeHouseholderInPlace(tau_v, diagonal[k]); if(k+10) A.coeffRef(k-1,k) = tau_u_prev; tau_u_prev = tau_u; } else A.coeffRef(k-1,k) = tau_u_prev; A.coeffRef(k,k) = tau_v; } if(bsbs && brows>bs) { SubMatType A22( A.bottomRightCorner(brows-bs,bcols-bs) ); SubMatType A21( A.block(bs,0, brows-bs,bs) ); SubMatType A12( A.block(0,bs, bs, bcols-bs) ); Scalar tmp = A12(bs-1,0); A12(bs-1,0) = 1; A22.noalias() -= A21 * Y.topLeftCorner(bcols,bs).bottomRows(bcols-bs).adjoint(); A22.noalias() -= X.topLeftCorner(brows,bs).bottomRows(brows-bs) * A12; A12(bs-1,0) = tmp; } } /** \internal * * Implementation of a block-bidiagonal reduction. * It is based on the following paper: * The Design of a Parallel Dense Linear Algebra Software Library: Reduction to Hessenberg, Tridiagonal, and Bidiagonal Form. * by Jaeyoung Choi, Jack J. Dongarra, David W. Walker. (1995) * section 3.3 */ template void upperbidiagonalization_inplace_blocked(MatrixType& A, BidiagType& bidiagonal, typename MatrixType::Index maxBlockSize=32, typename MatrixType::Scalar* tempData = 0) { typedef typename MatrixType::Index Index; typedef typename MatrixType::Scalar Scalar; typedef Block BlockType; Index rows = A.rows(); Index cols = A.cols(); Index size = (std::min)(rows, cols); typedef Matrix TempType; TempType tempVector; if(tempData==0) { tempVector.resize(cols); tempData = tempVector.data(); } Matrix X(rows,maxBlockSize), Y(cols, maxBlockSize); X.setOnes(); Y.setOnes(); Index blockSize = (std::min)(maxBlockSize,size); Index k = 0; for (k = 0; k < size; k += blockSize) { Index bs = (std::min)(size-k,blockSize); // actual size of the block Index tcols = cols - k - bs; // trailing columns Index trows = rows - k - bs; // trailing rows Index brows = rows - k; // rows of the block Index bcols = cols - k; // columns of the block // partition the matrix A: // // | A00 A01 A02 | // | | // A = | A10 A11 A12 | // | | // | A20 A21 A22 | // // where A11 is a bs x bs diagonal block, // and performs the bidiagonalization of A11, A21, A12, without updating A22. // // A22 will be updated in a second stage using matrices X and Y and level 3 operations: // A22 -= V*Y^T - X*U^T // where V and U contains the left and right Householder vectors // // Finally, the algorithm continue on the updated A22. // Let: // | A11 A12 | // B = | | // | A21 A22 | BlockType B = A.block(k,k,brows,bcols); upperbidiagonalization_inplace_helper(B, &(bidiagonal.template diagonal<0>().coeffRef(k)), &(bidiagonal.template diagonal<1>().coeffRef(k)), bs, X.topLeftCorner(brows,bs), Y.topLeftCorner(bcols,bs) ); } } // Standard upper bidiagonalization without fancy optimizations // This version should be faster for small matrix size template void upperbidiagonalization_inplace_unblocked(MatrixType& mat, BidiagType& bidiagonal, typename MatrixType::Scalar* tempData = 0) { typedef typename MatrixType::Index Index; typedef typename MatrixType::Scalar Scalar; Index rows = mat.rows(); Index cols = mat.cols(); Index size = (std::min)(rows, cols); typedef Matrix TempType; TempType tempVector; if(tempData==0) { tempVector.resize(rows); tempData = tempVector.data(); } for (Index k = 0; /* breaks at k==cols-1 below */ ; ++k) { Index remainingRows = rows - k; Index remainingCols = cols - k - 1; // construct left householder transform in-place in A mat.col(k).tail(remainingRows) .makeHouseholderInPlace(mat.coeffRef(k,k), bidiagonal.template diagonal<0>().coeffRef(k)); // apply householder transform to remaining part of A on the left mat.bottomRightCorner(remainingRows, remainingCols) .applyHouseholderOnTheLeft(mat.col(k).tail(remainingRows-1), mat.coeff(k,k), tempData); if(k == cols-1) break; // construct right householder transform in-place in mat mat.row(k).tail(remainingCols) .makeHouseholderInPlace(mat.coeffRef(k,k+1), bidiagonal.template diagonal<1>().coeffRef(k)); // apply householder transform to remaining part of mat on the left mat.bottomRightCorner(remainingRows-1, remainingCols) .applyHouseholderOnTheRight(mat.row(k).tail(remainingCols-1).transpose(), mat.coeff(k,k+1), tempData); } } template UpperBidiagonalization<_MatrixType>& UpperBidiagonalization<_MatrixType>::computeUnblocked(const _MatrixType& matrix) { Index rows = matrix.rows(); Index cols = matrix.cols(); eigen_assert(rows >= cols && "UpperBidiagonalization is only for Arices satisfying rows>=cols."); m_householder = matrix; ColVectorType temp(rows); upperbidiagonalization_inplace_unblocked(m_householder, m_bidiagonal, temp.data()); MatrixType A = matrix; BidiagonalType B(cols,cols); upperbidiagonalization_inplace_blocked(A, B, 8); std::cout << (m_householder-A)/*.maxCoeff()*/ << "\n\n"; std::cout << m_householder << "\n\n" << m_bidiagonal.template diagonal<0>() << "\n" << m_bidiagonal.template diagonal<1>() << "\n\n"; std::cout << A << "\n\n" << B.template diagonal<0>() << "\n" << B.template diagonal<1>() << "\n\n"; m_isInitialized = true; return *this; } template UpperBidiagonalization<_MatrixType>& UpperBidiagonalization<_MatrixType>::compute(const _MatrixType& matrix) { Index rows = matrix.rows(); Index cols = matrix.cols(); eigen_assert(rows >= cols && "UpperBidiagonalization is only for Arices satisfying rows>=cols."); m_householder = matrix; upperbidiagonalization_inplace_blocked(m_householder, m_bidiagonal); m_isInitialized = true; return *this; } #if 0 /** \return the Householder QR decomposition of \c *this. * * \sa class Bidiagonalization */ template const UpperBidiagonalization::PlainObject> MatrixBase::bidiagonalization() const { return UpperBidiagonalization(eval()); } #endif } // end namespace internal } // end namespace Eigen #endif // EIGEN_BIDIAGONALIZATION_H