aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/IterativeLinearSolvers10
-rw-r--r--Eigen/src/IterativeLinearSolvers/BasicPreconditioners.h70
-rw-r--r--Eigen/src/IterativeLinearSolvers/ConjugateGradient.h44
-rw-r--r--Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h213
-rw-r--r--test/CMakeLists.txt1
-rw-r--r--test/lscg.cpp29
-rw-r--r--test/sparse_solver.h58
7 files changed, 392 insertions, 33 deletions
diff --git a/Eigen/IterativeLinearSolvers b/Eigen/IterativeLinearSolvers
index c06668bd2..0594feb41 100644
--- a/Eigen/IterativeLinearSolvers
+++ b/Eigen/IterativeLinearSolvers
@@ -12,24 +12,26 @@
* This module currently provides iterative methods to solve problems of the form \c A \c x = \c b, where \c A is a squared matrix, usually very large and sparse.
* Those solvers are accessible via the following classes:
* - ConjugateGradient for selfadjoint (hermitian) matrices,
+ * - LSCG for rectangular least-square problems,
* - BiCGSTAB for general square matrices.
*
* These iterative solvers are associated with some preconditioners:
* - IdentityPreconditioner - not really useful
* - DiagonalPreconditioner - also called JAcobi preconditioner, work very well on diagonal dominant matrices.
- * - IncompleteILUT - incomplete LU factorization with dual thresholding
+ * - IncompleteLUT - incomplete LU factorization with dual thresholding
*
* Such problems can also be solved using the direct sparse decomposition modules: SparseCholesky, CholmodSupport, UmfPackSupport, SuperLUSupport.
*
- * \code
- * #include <Eigen/IterativeLinearSolvers>
- * \endcode
+ \code
+ #include <Eigen/IterativeLinearSolvers>
+ \endcode
*/
#include "src/IterativeLinearSolvers/SolveWithGuess.h"
#include "src/IterativeLinearSolvers/IterativeSolverBase.h"
#include "src/IterativeLinearSolvers/BasicPreconditioners.h"
#include "src/IterativeLinearSolvers/ConjugateGradient.h"
+#include "src/IterativeLinearSolvers/LeastSquareConjugateGradient.h"
#include "src/IterativeLinearSolvers/BiCGSTAB.h"
#include "src/IterativeLinearSolvers/IncompleteLUT.h"
diff --git a/Eigen/src/IterativeLinearSolvers/BasicPreconditioners.h b/Eigen/src/IterativeLinearSolvers/BasicPreconditioners.h
index a09f81225..6da423cf6 100644
--- a/Eigen/src/IterativeLinearSolvers/BasicPreconditioners.h
+++ b/Eigen/src/IterativeLinearSolvers/BasicPreconditioners.h
@@ -17,9 +17,9 @@ namespace Eigen {
*
* This class allows to approximately solve for A.x = b problems assuming A is a diagonal matrix.
* In other words, this preconditioner neglects all off diagonal entries and, in Eigen's language, solves for:
- * \code
- * A.diagonal().asDiagonal() . x = b
- * \endcode
+ \code
+ A.diagonal().asDiagonal() . x = b
+ \endcode
*
* \tparam _Scalar the type of the scalar.
*
@@ -28,6 +28,7 @@ namespace Eigen {
*
* \note A variant that has yet to be implemented would attempt to preserve the norm of each column.
*
+ * \sa class LeastSquareDiagonalPreconditioner, class ConjugateGradient
*/
template <typename _Scalar>
class DiagonalPreconditioner
@@ -100,6 +101,69 @@ class DiagonalPreconditioner
bool m_isInitialized;
};
+/** \ingroup IterativeLinearSolvers_Module
+ * \brief Jacobi preconditioner for LSCG
+ *
+ * This class allows to approximately solve for A' A x = A' b problems assuming A' A is a diagonal matrix.
+ * In other words, this preconditioner neglects all off diagonal entries and, in Eigen's language, solves for:
+ \code
+ (A.adjoint() * A).diagonal().asDiagonal() * x = b
+ \endcode
+ *
+ * \tparam _Scalar the type of the scalar.
+ *
+ * The diagonal entries are pre-inverted and stored into a dense vector.
+ *
+ * \sa class LSCG, class DiagonalPreconditioner
+ */
+template <typename _Scalar>
+class LeastSquareDiagonalPreconditioner : public DiagonalPreconditioner<_Scalar>
+{
+ typedef _Scalar Scalar;
+ typedef typename NumTraits<Scalar>::Real RealScalar;
+ typedef DiagonalPreconditioner<_Scalar> Base;
+ using Base::m_invdiag;
+ public:
+
+ LeastSquareDiagonalPreconditioner() : Base() {}
+
+ template<typename MatType>
+ explicit LeastSquareDiagonalPreconditioner(const MatType& mat) : Base()
+ {
+ compute(mat);
+ }
+
+ template<typename MatType>
+ LeastSquareDiagonalPreconditioner& analyzePattern(const MatType& )
+ {
+ return *this;
+ }
+
+ template<typename MatType>
+ LeastSquareDiagonalPreconditioner& factorize(const MatType& mat)
+ {
+ // Compute the inverse squared-norm of each column of mat
+ m_invdiag.resize(mat.cols());
+ for(Index j=0; j<mat.outerSize(); ++j)
+ {
+ RealScalar sum = mat.innerVector(j).squaredNorm();
+ if(sum>0)
+ m_invdiag(j) = RealScalar(1)/sum;
+ else
+ m_invdiag(j) = RealScalar(1);
+ }
+ Base::m_isInitialized = true;
+ return *this;
+ }
+
+ template<typename MatType>
+ LeastSquareDiagonalPreconditioner& compute(const MatType& mat)
+ {
+ return factorize(mat);
+ }
+
+ protected:
+};
/** \ingroup IterativeLinearSolvers_Module
* \brief A naive preconditioner which approximates any matrix as the identity matrix
diff --git a/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h b/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
index a799c3ef5..10cd94783 100644
--- a/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
+++ b/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
@@ -60,29 +60,29 @@ void conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x,
}
VectorType p(n);
- p = precond.solve(residual); //initial search direction
+ p = precond.solve(residual); // initial search direction
VectorType z(n), tmp(n);
RealScalar absNew = numext::real(residual.dot(p)); // the square of the absolute value of r scaled by invM
Index i = 0;
while(i < maxIters)
{
- tmp.noalias() = mat * p; // the bottleneck of the algorithm
+ tmp.noalias() = mat * p; // the bottleneck of the algorithm
- Scalar alpha = absNew / p.dot(tmp); // the amount we travel on dir
- x += alpha * p; // update solution
- residual -= alpha * tmp; // update residue
+ Scalar alpha = absNew / p.dot(tmp); // the amount we travel on dir
+ x += alpha * p; // update solution
+ residual -= alpha * tmp; // update residual
residualNorm2 = residual.squaredNorm();
if(residualNorm2 < threshold)
break;
- z = precond.solve(residual); // approximately solve for "A z = residual"
+ z = precond.solve(residual); // approximately solve for "A z = residual"
RealScalar absOld = absNew;
absNew = numext::real(residual.dot(z)); // update the absolute value of r
- RealScalar beta = absNew / absOld; // calculate the Gram-Schmidt value used to create the new search direction
- p = z + beta * p; // update search direction
+ RealScalar beta = absNew / absOld; // calculate the Gram-Schmidt value used to create the new search direction
+ p = z + beta * p; // update search direction
i++;
}
tol_error = sqrt(residualNorm2 / rhsNorm2);
@@ -122,24 +122,24 @@ struct traits<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
* and NumTraits<Scalar>::epsilon() for the tolerance.
*
* This class can be used as the direct solver classes. Here is a typical usage example:
- * \code
- * int n = 10000;
- * VectorXd x(n), b(n);
- * SparseMatrix<double> A(n,n);
- * // fill A and b
- * ConjugateGradient<SparseMatrix<double> > cg;
- * cg.compute(A);
- * x = cg.solve(b);
- * std::cout << "#iterations: " << cg.iterations() << std::endl;
- * std::cout << "estimated error: " << cg.error() << std::endl;
- * // update b, and solve again
- * x = cg.solve(b);
- * \endcode
+ \code
+ int n = 10000;
+ VectorXd x(n), b(n);
+ SparseMatrix<double> A(n,n);
+ // fill A and b
+ ConjugateGradient<SparseMatrix<double> > cg;
+ cg.compute(A);
+ x = cg.solve(b);
+ std::cout << "#iterations: " << cg.iterations() << std::endl;
+ std::cout << "estimated error: " << cg.error() << std::endl;
+ // update b, and solve again
+ x = cg.solve(b);
+ \endcode
*
* By default the iterations start with x=0 as an initial guess of the solution.
* One can control the start using the solveWithGuess() method.
*
- * \sa class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
+ * \sa class LSCG, class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
*/
template< typename _MatrixType, int _UpLo, typename _Preconditioner>
class ConjugateGradient : public IterativeSolverBase<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
diff --git a/Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h b/Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h
new file mode 100644
index 000000000..beaf5c307
--- /dev/null
+++ b/Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h
@@ -0,0 +1,213 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2015 Gael Guennebaud <gael.guennebaud@inria.fr>
+//
+// 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_LEAST_SQUARE_CONJUGATE_GRADIENT_H
+#define EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
+
+namespace Eigen {
+
+namespace internal {
+
+/** \internal Low-level conjugate gradient algorithm for least-square problems
+ * \param mat The matrix A
+ * \param rhs The right hand side vector b
+ * \param x On input and initial solution, on output the computed solution.
+ * \param precond A preconditioner being able to efficiently solve for an
+ * approximation of A'Ax=b (regardless of b)
+ * \param iters On input the max number of iteration, on output the number of performed iterations.
+ * \param tol_error On input the tolerance error, on output an estimation of the relative error.
+ */
+template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
+EIGEN_DONT_INLINE
+void least_square_conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x,
+ const Preconditioner& precond, Index& iters,
+ typename Dest::RealScalar& tol_error)
+{
+ using std::sqrt;
+ using std::abs;
+ typedef typename Dest::RealScalar RealScalar;
+ typedef typename Dest::Scalar Scalar;
+ typedef Matrix<Scalar,Dynamic,1> VectorType;
+
+ RealScalar tol = tol_error;
+ Index maxIters = iters;
+
+ Index m = mat.rows(), n = mat.cols();
+
+ VectorType residual = rhs - mat * x;
+ VectorType normal_residual = mat.adjoint() * residual;
+
+ RealScalar rhsNorm2 = (mat.adjoint()*rhs).squaredNorm();
+ if(rhsNorm2 == 0)
+ {
+ x.setZero();
+ iters = 0;
+ tol_error = 0;
+ return;
+ }
+ RealScalar threshold = tol*tol*rhsNorm2;
+ RealScalar residualNorm2 = normal_residual.squaredNorm();
+ if (residualNorm2 < threshold)
+ {
+ iters = 0;
+ tol_error = sqrt(residualNorm2 / rhsNorm2);
+ return;
+ }
+
+ VectorType p(n);
+ p = precond.solve(normal_residual); // initial search direction
+
+ VectorType z(n), tmp(m);
+ RealScalar absNew = numext::real(normal_residual.dot(p)); // the square of the absolute value of r scaled by invM
+ Index i = 0;
+ while(i < maxIters)
+ {
+ tmp.noalias() = mat * p;
+
+ Scalar alpha = absNew / tmp.squaredNorm(); // the amount we travel on dir
+ x += alpha * p; // update solution
+ residual -= alpha * tmp; // update residual
+ normal_residual = mat.adjoint() * residual; // update residual of the normal equation
+
+ residualNorm2 = normal_residual.squaredNorm();
+ if(residualNorm2 < threshold)
+ break;
+
+ z = precond.solve(normal_residual); // approximately solve for "A'A z = normal_residual"
+
+ RealScalar absOld = absNew;
+ absNew = numext::real(normal_residual.dot(z)); // update the absolute value of r
+ RealScalar beta = absNew / absOld; // calculate the Gram-Schmidt value used to create the new search direction
+ p = z + beta * p; // update search direction
+ i++;
+ }
+ tol_error = sqrt(residualNorm2 / rhsNorm2);
+ iters = i;
+}
+
+}
+
+template< typename _MatrixType,
+ typename _Preconditioner = LeastSquareDiagonalPreconditioner<typename _MatrixType::Scalar> >
+class LSCG;
+
+namespace internal {
+
+template< typename _MatrixType, typename _Preconditioner>
+struct traits<LSCG<_MatrixType,_Preconditioner> >
+{
+ typedef _MatrixType MatrixType;
+ typedef _Preconditioner Preconditioner;
+};
+
+}
+
+/** \ingroup IterativeLinearSolvers_Module
+ * \brief A conjugate gradient solver for sparse (or dense) least-square problems
+ *
+ * This class allows to solve for A x = b linear problems using an iterative conjugate gradient algorithm.
+ * The matrix A can be non symmetric and rectangular, but the matrix A' A should be positive-definite to guaranty stability.
+ * Otherwise, the SparseLU or SparseQR classes might be preferable.
+ * The matrix A and the vectors x and b can be either dense or sparse.
+ *
+ * \tparam _MatrixType the type of the matrix A, can be a dense or a sparse matrix.
+ * \tparam _Preconditioner the type of the preconditioner. Default is LeastSquareDiagonalPreconditioner
+ *
+ * The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
+ * and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
+ * and NumTraits<Scalar>::epsilon() for the tolerance.
+ *
+ * This class can be used as the direct solver classes. Here is a typical usage example:
+ \code
+ int m=1000000, n = 10000;
+ VectorXd x(n), b(m);
+ SparseMatrix<double> A(m,n);
+ // fill A and b
+ LSCG<SparseMatrix<double> > lscg;
+ lscg.compute(A);
+ x = lscg.solve(b);
+ std::cout << "#iterations: " << lscg.iterations() << std::endl;
+ std::cout << "estimated error: " << lscg.error() << std::endl;
+ // update b, and solve again
+ x = lscg.solve(b);
+ \endcode
+ *
+ * By default the iterations start with x=0 as an initial guess of the solution.
+ * One can control the start using the solveWithGuess() method.
+ *
+ * \sa class ConjugateGradient, SparseLU, SparseQR
+ */
+template< typename _MatrixType, typename _Preconditioner>
+class LSCG : public IterativeSolverBase<LSCG<_MatrixType,_Preconditioner> >
+{
+ typedef IterativeSolverBase<LSCG> Base;
+ using Base::mp_matrix;
+ using Base::m_error;
+ using Base::m_iterations;
+ using Base::m_info;
+ using Base::m_isInitialized;
+public:
+ typedef _MatrixType MatrixType;
+ typedef typename MatrixType::Scalar Scalar;
+ typedef typename MatrixType::RealScalar RealScalar;
+ typedef _Preconditioner Preconditioner;
+
+public:
+
+ /** Default constructor. */
+ LSCG() : Base() {}
+
+ /** Initialize the solver with matrix \a A for further \c Ax=b solving.
+ *
+ * This constructor is a shortcut for the default constructor followed
+ * by a call to compute().
+ *
+ * \warning this class stores a reference to the matrix A as well as some
+ * precomputed values that depend on it. Therefore, if \a A is changed
+ * this class becomes invalid. Call compute() to update it with the new
+ * matrix A, or modify a copy of A.
+ */
+ explicit LSCG(const MatrixType& A) : Base(A) {}
+
+ ~LSCG() {}
+
+ /** \internal */
+ template<typename Rhs,typename Dest>
+ void _solve_with_guess_impl(const Rhs& b, Dest& x) const
+ {
+ m_iterations = Base::maxIterations();
+ m_error = Base::m_tolerance;
+
+ for(Index j=0; j<b.cols(); ++j)
+ {
+ m_iterations = Base::maxIterations();
+ m_error = Base::m_tolerance;
+
+ typename Dest::ColXpr xj(x,j);
+ internal::least_square_conjugate_gradient(mp_matrix, b.col(j), xj, Base::m_preconditioner, m_iterations, m_error);
+ }
+
+ m_isInitialized = true;
+ m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
+ }
+
+ /** \internal */
+ using Base::_solve_impl;
+ template<typename Rhs,typename Dest>
+ void _solve_impl(const MatrixBase<Rhs>& b, Dest& x) const
+ {
+ x.setZero();
+ _solve_with_guess_impl(b.derived(),x);
+ }
+
+};
+
+} // end namespace Eigen
+
+#endif // EIGEN_LEAST_SQUARE_CONJUGATE_GRADIENT_H
diff --git a/test/CMakeLists.txt b/test/CMakeLists.txt
index 168749634..1712b8718 100644
--- a/test/CMakeLists.txt
+++ b/test/CMakeLists.txt
@@ -234,6 +234,7 @@ ei_add_test(sparse_permutations)
ei_add_test(simplicial_cholesky)
ei_add_test(conjugate_gradient)
ei_add_test(bicgstab)
+ei_add_test(lscg)
ei_add_test(sparselu)
ei_add_test(sparseqr)
ei_add_test(umeyama)
diff --git a/test/lscg.cpp b/test/lscg.cpp
new file mode 100644
index 000000000..599ed5619
--- /dev/null
+++ b/test/lscg.cpp
@@ -0,0 +1,29 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2011 Gael Guennebaud <g.gael@free.fr>
+//
+// 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/.
+
+#include "sparse_solver.h"
+#include <Eigen/IterativeLinearSolvers>
+
+template<typename T> void test_lscg_T()
+{
+ LSCG<SparseMatrix<T> > lscg_colmajor_diag;
+ LSCG<SparseMatrix<T>, IdentityPreconditioner> lscg_colmajor_I;
+
+ CALL_SUBTEST( check_sparse_square_solving(lscg_colmajor_diag) );
+ CALL_SUBTEST( check_sparse_square_solving(lscg_colmajor_I) );
+
+ CALL_SUBTEST( check_sparse_leastsquare_solving(lscg_colmajor_diag) );
+ CALL_SUBTEST( check_sparse_leastsquare_solving(lscg_colmajor_I) );
+}
+
+void test_lscg()
+{
+ CALL_SUBTEST_1(test_lscg_T<double>());
+ CALL_SUBTEST_2(test_lscg_T<std::complex<double> >());
+}
diff --git a/test/sparse_solver.h b/test/sparse_solver.h
index f0a4691e3..f266e2c9a 100644
--- a/test/sparse_solver.h
+++ b/test/sparse_solver.h
@@ -17,9 +17,9 @@ void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A,
typedef typename Mat::Scalar Scalar;
typedef typename Mat::StorageIndex StorageIndex;
- DenseRhs refX = dA.lu().solve(db);
+ DenseRhs refX = dA.householderQr().solve(db);
{
- Rhs x(b.rows(), b.cols());
+ Rhs x(A.cols(), b.cols());
Rhs oldb = b;
solver.compute(A);
@@ -94,7 +94,7 @@ void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A,
// test dense Block as the result and rhs:
{
- DenseRhs x(db.rows(), db.cols());
+ DenseRhs x(refX.rows(), refX.cols());
DenseRhs oldb(db);
x.setZero();
x.block(0,0,x.rows(),x.cols()) = solver.solve(db.block(0,0,db.rows(),db.cols()));
@@ -119,7 +119,7 @@ void check_sparse_solving_real_cases(Solver& solver, const typename Solver::Matr
typedef typename Mat::Scalar Scalar;
typedef typename Mat::RealScalar RealScalar;
- Rhs x(b.rows(), b.cols());
+ Rhs x(A.cols(), b.cols());
solver.compute(A);
if (solver.info() != Success)
@@ -410,3 +410,53 @@ template<typename Solver> void check_sparse_square_abs_determinant(Solver& solve
}
}
+template<typename Solver, typename DenseMat>
+void generate_sparse_leastsquare_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300, int options = ForceNonZeroDiag)
+{
+ typedef typename Solver::MatrixType Mat;
+ typedef typename Mat::Scalar Scalar;
+
+ int rows = internal::random<int>(1,maxSize);
+ int cols = internal::random<int>(1,rows);
+ double density = (std::max)(8./(rows*cols), 0.01);
+
+ A.resize(rows,cols);
+ dA.resize(rows,cols);
+
+ initSparse<Scalar>(density, dA, A, options);
+}
+
+template<typename Solver> void check_sparse_leastsquare_solving(Solver& solver)
+{
+ typedef typename Solver::MatrixType Mat;
+ typedef typename Mat::Scalar Scalar;
+ typedef SparseMatrix<Scalar,ColMajor> SpMat;
+ typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
+ typedef Matrix<Scalar,Dynamic,1> DenseVector;
+
+ int rhsCols = internal::random<int>(1,16);
+
+ Mat A;
+ DenseMatrix dA;
+ for (int i = 0; i < g_repeat; i++) {
+ generate_sparse_leastsquare_problem(solver, A, dA);
+
+ A.makeCompressed();
+ DenseVector b = DenseVector::Random(A.rows());
+ DenseMatrix dB(A.rows(),rhsCols);
+ SpMat B(A.rows(),rhsCols);
+ double density = (std::max)(8./(A.rows()*rhsCols), 0.1);
+ initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
+ B.makeCompressed();
+ check_sparse_solving(solver, A, b, dA, b);
+ check_sparse_solving(solver, A, dB, dA, dB);
+ check_sparse_solving(solver, A, B, dA, dB);
+
+ // check only once
+ if(i==0)
+ {
+ b = DenseVector::Zero(A.rows());
+ check_sparse_solving(solver, A, b, dA, b);
+ }
+ }
+}