blob: f7ce471349a2b39c57d2f0ad8f425c512f99cc79 [file] [log] [blame]
Brian Silverman72890c22015-09-19 14:37:37 -04001// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
Austin Schuh189376f2018-12-20 22:11:15 +11004// Copyright (C) 2011-2014 Gael Guennebaud <gael.guennebaud@inria.fr>
Brian Silverman72890c22015-09-19 14:37:37 -04005//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_CONJUGATE_GRADIENT_H
11#define EIGEN_CONJUGATE_GRADIENT_H
12
13namespace Eigen {
14
15namespace internal {
16
17/** \internal Low-level conjugate gradient algorithm
18 * \param mat The matrix A
19 * \param rhs The right hand side vector b
20 * \param x On input and initial solution, on output the computed solution.
21 * \param precond A preconditioner being able to efficiently solve for an
22 * approximation of Ax=b (regardless of b)
23 * \param iters On input the max number of iteration, on output the number of performed iterations.
24 * \param tol_error On input the tolerance error, on output an estimation of the relative error.
25 */
26template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
27EIGEN_DONT_INLINE
28void conjugate_gradient(const MatrixType& mat, const Rhs& rhs, Dest& x,
Austin Schuh189376f2018-12-20 22:11:15 +110029 const Preconditioner& precond, Index& iters,
Brian Silverman72890c22015-09-19 14:37:37 -040030 typename Dest::RealScalar& tol_error)
31{
32 using std::sqrt;
33 using std::abs;
34 typedef typename Dest::RealScalar RealScalar;
35 typedef typename Dest::Scalar Scalar;
36 typedef Matrix<Scalar,Dynamic,1> VectorType;
37
38 RealScalar tol = tol_error;
Austin Schuh189376f2018-12-20 22:11:15 +110039 Index maxIters = iters;
Brian Silverman72890c22015-09-19 14:37:37 -040040
Austin Schuh189376f2018-12-20 22:11:15 +110041 Index n = mat.cols();
Brian Silverman72890c22015-09-19 14:37:37 -040042
43 VectorType residual = rhs - mat * x; //initial residual
44
45 RealScalar rhsNorm2 = rhs.squaredNorm();
46 if(rhsNorm2 == 0)
47 {
48 x.setZero();
49 iters = 0;
50 tol_error = 0;
51 return;
52 }
Austin Schuh189376f2018-12-20 22:11:15 +110053 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
54 RealScalar threshold = numext::maxi(tol*tol*rhsNorm2,considerAsZero);
Brian Silverman72890c22015-09-19 14:37:37 -040055 RealScalar residualNorm2 = residual.squaredNorm();
56 if (residualNorm2 < threshold)
57 {
58 iters = 0;
59 tol_error = sqrt(residualNorm2 / rhsNorm2);
60 return;
61 }
Austin Schuh189376f2018-12-20 22:11:15 +110062
Brian Silverman72890c22015-09-19 14:37:37 -040063 VectorType p(n);
Austin Schuh189376f2018-12-20 22:11:15 +110064 p = precond.solve(residual); // initial search direction
Brian Silverman72890c22015-09-19 14:37:37 -040065
66 VectorType z(n), tmp(n);
67 RealScalar absNew = numext::real(residual.dot(p)); // the square of the absolute value of r scaled by invM
Austin Schuh189376f2018-12-20 22:11:15 +110068 Index i = 0;
Brian Silverman72890c22015-09-19 14:37:37 -040069 while(i < maxIters)
70 {
Austin Schuh189376f2018-12-20 22:11:15 +110071 tmp.noalias() = mat * p; // the bottleneck of the algorithm
Brian Silverman72890c22015-09-19 14:37:37 -040072
Austin Schuh189376f2018-12-20 22:11:15 +110073 Scalar alpha = absNew / p.dot(tmp); // the amount we travel on dir
74 x += alpha * p; // update solution
75 residual -= alpha * tmp; // update residual
Brian Silverman72890c22015-09-19 14:37:37 -040076
77 residualNorm2 = residual.squaredNorm();
78 if(residualNorm2 < threshold)
79 break;
80
Austin Schuh189376f2018-12-20 22:11:15 +110081 z = precond.solve(residual); // approximately solve for "A z = residual"
Brian Silverman72890c22015-09-19 14:37:37 -040082
83 RealScalar absOld = absNew;
84 absNew = numext::real(residual.dot(z)); // update the absolute value of r
Austin Schuh189376f2018-12-20 22:11:15 +110085 RealScalar beta = absNew / absOld; // calculate the Gram-Schmidt value used to create the new search direction
86 p = z + beta * p; // update search direction
Brian Silverman72890c22015-09-19 14:37:37 -040087 i++;
88 }
89 tol_error = sqrt(residualNorm2 / rhsNorm2);
90 iters = i;
91}
92
93}
94
95template< typename _MatrixType, int _UpLo=Lower,
96 typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar> >
97class ConjugateGradient;
98
99namespace internal {
100
101template< typename _MatrixType, int _UpLo, typename _Preconditioner>
102struct traits<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
103{
104 typedef _MatrixType MatrixType;
105 typedef _Preconditioner Preconditioner;
106};
107
108}
109
110/** \ingroup IterativeLinearSolvers_Module
Austin Schuh189376f2018-12-20 22:11:15 +1100111 * \brief A conjugate gradient solver for sparse (or dense) self-adjoint problems
Brian Silverman72890c22015-09-19 14:37:37 -0400112 *
Austin Schuh189376f2018-12-20 22:11:15 +1100113 * This class allows to solve for A.x = b linear problems using an iterative conjugate gradient algorithm.
114 * The matrix A must be selfadjoint. The matrix A and the vectors x and b can be either dense or sparse.
Brian Silverman72890c22015-09-19 14:37:37 -0400115 *
116 * \tparam _MatrixType the type of the matrix A, can be a dense or a sparse matrix.
117 * \tparam _UpLo the triangular part that will be used for the computations. It can be Lower,
Austin Schuh189376f2018-12-20 22:11:15 +1100118 * \c Upper, or \c Lower|Upper in which the full matrix entries will be considered.
119 * Default is \c Lower, best performance is \c Lower|Upper.
Brian Silverman72890c22015-09-19 14:37:37 -0400120 * \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
121 *
Austin Schuh189376f2018-12-20 22:11:15 +1100122 * \implsparsesolverconcept
123 *
Brian Silverman72890c22015-09-19 14:37:37 -0400124 * The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
125 * and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations
126 * and NumTraits<Scalar>::epsilon() for the tolerance.
127 *
Austin Schuh189376f2018-12-20 22:11:15 +1100128 * The tolerance corresponds to the relative residual error: |Ax-b|/|b|
129 *
130 * \b Performance: Even though the default value of \c _UpLo is \c Lower, significantly higher performance is
131 * achieved when using a complete matrix and \b Lower|Upper as the \a _UpLo template parameter. Moreover, in this
132 * case multi-threading can be exploited if the user code is compiled with OpenMP enabled.
133 * See \ref TopicMultiThreading for details.
134 *
Brian Silverman72890c22015-09-19 14:37:37 -0400135 * This class can be used as the direct solver classes. Here is a typical usage example:
Austin Schuh189376f2018-12-20 22:11:15 +1100136 \code
137 int n = 10000;
138 VectorXd x(n), b(n);
139 SparseMatrix<double> A(n,n);
140 // fill A and b
141 ConjugateGradient<SparseMatrix<double>, Lower|Upper> cg;
142 cg.compute(A);
143 x = cg.solve(b);
144 std::cout << "#iterations: " << cg.iterations() << std::endl;
145 std::cout << "estimated error: " << cg.error() << std::endl;
146 // update b, and solve again
147 x = cg.solve(b);
148 \endcode
Brian Silverman72890c22015-09-19 14:37:37 -0400149 *
150 * By default the iterations start with x=0 as an initial guess of the solution.
151 * One can control the start using the solveWithGuess() method.
152 *
Austin Schuh189376f2018-12-20 22:11:15 +1100153 * ConjugateGradient can also be used in a matrix-free context, see the following \link MatrixfreeSolverExample example \endlink.
154 *
155 * \sa class LeastSquaresConjugateGradient, class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
Brian Silverman72890c22015-09-19 14:37:37 -0400156 */
157template< typename _MatrixType, int _UpLo, typename _Preconditioner>
158class ConjugateGradient : public IterativeSolverBase<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
159{
160 typedef IterativeSolverBase<ConjugateGradient> Base;
Austin Schuh189376f2018-12-20 22:11:15 +1100161 using Base::matrix;
Brian Silverman72890c22015-09-19 14:37:37 -0400162 using Base::m_error;
163 using Base::m_iterations;
164 using Base::m_info;
165 using Base::m_isInitialized;
166public:
167 typedef _MatrixType MatrixType;
168 typedef typename MatrixType::Scalar Scalar;
Brian Silverman72890c22015-09-19 14:37:37 -0400169 typedef typename MatrixType::RealScalar RealScalar;
170 typedef _Preconditioner Preconditioner;
171
172 enum {
173 UpLo = _UpLo
174 };
175
176public:
177
178 /** Default constructor. */
179 ConjugateGradient() : Base() {}
180
181 /** Initialize the solver with matrix \a A for further \c Ax=b solving.
182 *
183 * This constructor is a shortcut for the default constructor followed
184 * by a call to compute().
185 *
186 * \warning this class stores a reference to the matrix A as well as some
187 * precomputed values that depend on it. Therefore, if \a A is changed
188 * this class becomes invalid. Call compute() to update it with the new
189 * matrix A, or modify a copy of A.
190 */
Austin Schuh189376f2018-12-20 22:11:15 +1100191 template<typename MatrixDerived>
192 explicit ConjugateGradient(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
Brian Silverman72890c22015-09-19 14:37:37 -0400193
194 ~ConjugateGradient() {}
Brian Silverman72890c22015-09-19 14:37:37 -0400195
196 /** \internal */
197 template<typename Rhs,typename Dest>
Austin Schuh189376f2018-12-20 22:11:15 +1100198 void _solve_with_guess_impl(const Rhs& b, Dest& x) const
Brian Silverman72890c22015-09-19 14:37:37 -0400199 {
Austin Schuh189376f2018-12-20 22:11:15 +1100200 typedef typename Base::MatrixWrapper MatrixWrapper;
201 typedef typename Base::ActualMatrixType ActualMatrixType;
202 enum {
203 TransposeInput = (!MatrixWrapper::MatrixFree)
204 && (UpLo==(Lower|Upper))
205 && (!MatrixType::IsRowMajor)
206 && (!NumTraits<Scalar>::IsComplex)
207 };
208 typedef typename internal::conditional<TransposeInput,Transpose<const ActualMatrixType>, ActualMatrixType const&>::type RowMajorWrapper;
209 EIGEN_STATIC_ASSERT(EIGEN_IMPLIES(MatrixWrapper::MatrixFree,UpLo==(Lower|Upper)),MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY);
Brian Silverman72890c22015-09-19 14:37:37 -0400210 typedef typename internal::conditional<UpLo==(Lower|Upper),
Austin Schuh189376f2018-12-20 22:11:15 +1100211 RowMajorWrapper,
212 typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type
213 >::type SelfAdjointWrapper;
Brian Silverman72890c22015-09-19 14:37:37 -0400214 m_iterations = Base::maxIterations();
215 m_error = Base::m_tolerance;
216
Austin Schuh189376f2018-12-20 22:11:15 +1100217 for(Index j=0; j<b.cols(); ++j)
Brian Silverman72890c22015-09-19 14:37:37 -0400218 {
219 m_iterations = Base::maxIterations();
220 m_error = Base::m_tolerance;
221
222 typename Dest::ColXpr xj(x,j);
Austin Schuh189376f2018-12-20 22:11:15 +1100223 RowMajorWrapper row_mat(matrix());
224 internal::conjugate_gradient(SelfAdjointWrapper(row_mat), b.col(j), xj, Base::m_preconditioner, m_iterations, m_error);
Brian Silverman72890c22015-09-19 14:37:37 -0400225 }
226
227 m_isInitialized = true;
228 m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
229 }
230
231 /** \internal */
Austin Schuh189376f2018-12-20 22:11:15 +1100232 using Base::_solve_impl;
Brian Silverman72890c22015-09-19 14:37:37 -0400233 template<typename Rhs,typename Dest>
Austin Schuh189376f2018-12-20 22:11:15 +1100234 void _solve_impl(const MatrixBase<Rhs>& b, Dest& x) const
Brian Silverman72890c22015-09-19 14:37:37 -0400235 {
236 x.setZero();
Austin Schuh189376f2018-12-20 22:11:15 +1100237 _solve_with_guess_impl(b.derived(),x);
Brian Silverman72890c22015-09-19 14:37:37 -0400238 }
239
240protected:
241
242};
243
Brian Silverman72890c22015-09-19 14:37:37 -0400244} // end namespace Eigen
245
246#endif // EIGEN_CONJUGATE_GRADIENT_H