$darkmode
Eigen-unsupported  5.0.1-dev
MINRES.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012 Giacomo Po <gpo@ucla.edu>
5 // Copyright (C) 2011-2014 Gael Guennebaud <gael.guennebaud@inria.fr>
6 // Copyright (C) 2018 David Hyde <dabh@stanford.edu>
7 //
8 // This Source Code Form is subject to the terms of the Mozilla
9 // Public License v. 2.0. If a copy of the MPL was not distributed
10 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
11 
12 #ifndef EIGEN_MINRES_H
13 #define EIGEN_MINRES_H
14 
15 // IWYU pragma: private
16 #include "./InternalHeaderCheck.h"
17 
18 namespace Eigen {
19 
20 namespace internal {
21 
31 template <typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
32 EIGEN_DONT_INLINE void minres(const MatrixType& mat, const Rhs& rhs, Dest& x, const Preconditioner& precond,
33  Index& iters, typename Dest::RealScalar& tol_error) {
34  using std::sqrt;
35  typedef typename Dest::RealScalar RealScalar;
36  typedef typename Dest::Scalar Scalar;
37  typedef Matrix<Scalar, Dynamic, 1> VectorType;
38 
39  // Check for zero rhs
40  const RealScalar rhsNorm2(rhs.squaredNorm());
41  if (rhsNorm2 == 0) {
42  x.setZero();
43  iters = 0;
44  tol_error = 0;
45  return;
46  }
47 
48  // initialize
49  const Index maxIters(iters); // initialize maxIters to iters
50  const Index N(mat.cols()); // the size of the matrix
51  const RealScalar threshold2(tol_error * tol_error * rhsNorm2); // convergence threshold (compared to residualNorm2)
52 
53  // Initialize preconditioned Lanczos
54  VectorType v_old(N); // will be initialized inside loop
55  VectorType v(VectorType::Zero(N)); // initialize v
56  VectorType v_new(rhs - mat * x); // initialize v_new
57  RealScalar residualNorm2(v_new.squaredNorm());
58  VectorType w(N); // will be initialized inside loop
59  VectorType w_new(precond.solve(v_new)); // initialize w_new
60  // RealScalar beta; // will be initialized inside loop
61  RealScalar beta_new2(v_new.dot(w_new));
62  eigen_assert(beta_new2 >= 0.0 && "PRECONDITIONER IS NOT POSITIVE DEFINITE");
63  RealScalar beta_new(sqrt(beta_new2));
64  const RealScalar beta_one(beta_new);
65  // Initialize other variables
66  RealScalar c(1.0); // the cosine of the Givens rotation
67  RealScalar c_old(1.0);
68  RealScalar s(0.0); // the sine of the Givens rotation
69  RealScalar s_old(0.0); // the sine of the Givens rotation
70  VectorType p_oold(N); // will be initialized in loop
71  VectorType p_old(VectorType::Zero(N)); // initialize p_old=0
72  VectorType p(p_old); // initialize p=0
73  RealScalar eta(1.0);
74 
75  iters = 0; // reset iters
76  while (iters < maxIters) {
77  // Preconditioned Lanczos
78  /* Note that there are 4 variants on the Lanczos algorithm. These are
79  * described in Paige, C. C. (1972). Computational variants of
80  * the Lanczos method for the eigenproblem. IMA Journal of Applied
81  * Mathematics, 10(3), 373-381. The current implementation corresponds
82  * to the case A(2,7) in the paper. It also corresponds to
83  * algorithm 6.14 in Y. Saad, Iterative Methods for Sparse Linear
84  * Systems, 2003 p.173. For the preconditioned version see
85  * A. Greenbaum, Iterative Methods for Solving Linear Systems, SIAM (1987).
86  */
87  const RealScalar beta(beta_new);
88  v_old = v; // update: at first time step, this makes v_old = 0 so value of beta doesn't matter
89  v_new /= beta_new; // overwrite v_new for next iteration
90  w_new /= beta_new; // overwrite w_new for next iteration
91  v = v_new; // update
92  w = w_new; // update
93  v_new.noalias() = mat * w - beta * v_old; // compute v_new
94  const RealScalar alpha = v_new.dot(w);
95  v_new -= alpha * v; // overwrite v_new
96  w_new = precond.solve(v_new); // overwrite w_new
97  beta_new2 = v_new.dot(w_new); // compute beta_new
98  eigen_assert(beta_new2 >= 0.0 && "PRECONDITIONER IS NOT POSITIVE DEFINITE");
99  beta_new = sqrt(beta_new2); // compute beta_new
100 
101  // Givens rotation
102  const RealScalar r2 = s * alpha + c * c_old * beta; // s, s_old, c and c_old are still from previous iteration
103  const RealScalar r3 = s_old * beta; // s, s_old, c and c_old are still from previous iteration
104  const RealScalar r1_hat = c * alpha - c_old * s * beta;
105  const RealScalar r1 = sqrt(std::pow(r1_hat, 2) + std::pow(beta_new, 2));
106  c_old = c; // store for next iteration
107  s_old = s; // store for next iteration
108  c = r1_hat / r1; // new cosine
109  s = beta_new / r1; // new sine
110 
111  // Update solution
112  p_oold = p_old;
113  p_old = p;
114  p.noalias() = (w - r2 * p_old - r3 * p_oold) / r1; // IS NOALIAS REQUIRED?
115  x += beta_one * c * eta * p;
116 
117  /* Update the squared residual. Note that this is the estimated residual.
118  The real residual |Ax-b|^2 may be slightly larger */
119  residualNorm2 *= s * s;
120 
121  if (residualNorm2 < threshold2) {
122  break;
123  }
124 
125  eta = -s * eta; // update eta
126  iters++; // increment iteration number (for output purposes)
127  }
128 
129  /* Compute error. Note that this is the estimated error. The real
130  error |Ax-b|/|b| may be slightly larger */
131  tol_error = std::sqrt(residualNorm2 / rhsNorm2);
132 }
133 
134 } // namespace internal
135 
136 template <typename MatrixType_, int UpLo_ = Lower, typename Preconditioner_ = IdentityPreconditioner>
137 class MINRES;
138 
139 namespace internal {
140 
141 template <typename MatrixType_, int UpLo_, typename Preconditioner_>
142 struct traits<MINRES<MatrixType_, UpLo_, Preconditioner_> > {
143  typedef MatrixType_ MatrixType;
144  typedef Preconditioner_ Preconditioner;
145 };
146 
147 } // namespace internal
148 
187 template <typename MatrixType_, int UpLo_, typename Preconditioner_>
188 class MINRES : public IterativeSolverBase<MINRES<MatrixType_, UpLo_, Preconditioner_> > {
189  typedef IterativeSolverBase<MINRES> Base;
190  using Base::m_error;
191  using Base::m_info;
192  using Base::m_isInitialized;
193  using Base::m_iterations;
194  using Base::matrix;
195 
196  public:
197  using Base::_solve_impl;
198  typedef MatrixType_ MatrixType;
199  typedef typename MatrixType::Scalar Scalar;
200  typedef typename MatrixType::RealScalar RealScalar;
201  typedef Preconditioner_ Preconditioner;
202 
203  enum { UpLo = UpLo_ };
204 
205  public:
207  MINRES() : Base() {}
208 
219  template <typename MatrixDerived>
220  explicit MINRES(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
221 
223  ~MINRES() {}
224 
226  template <typename Rhs, typename Dest>
227  void _solve_vector_with_guess_impl(const Rhs& b, Dest& x) const {
228  typedef typename Base::MatrixWrapper MatrixWrapper;
229  typedef typename Base::ActualMatrixType ActualMatrixType;
230  enum {
231  TransposeInput = (!MatrixWrapper::MatrixFree) && (UpLo == (Lower | Upper)) && (!MatrixType::IsRowMajor) &&
233  };
234  typedef std::conditional_t<TransposeInput, Transpose<const ActualMatrixType>, ActualMatrixType const&>
235  RowMajorWrapper;
236  EIGEN_STATIC_ASSERT(internal::check_implication(MatrixWrapper::MatrixFree, UpLo == (Lower | Upper)),
237  MATRIX_FREE_CONJUGATE_GRADIENT_IS_COMPATIBLE_WITH_UPPER_UNION_LOWER_MODE_ONLY);
238  typedef std::conditional_t<UpLo == (Lower | Upper), RowMajorWrapper,
239  typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type>
240  SelfAdjointWrapper;
241 
242  m_iterations = Base::maxIterations();
243  m_error = Base::m_tolerance;
244  RowMajorWrapper row_mat(matrix());
245  internal::minres(SelfAdjointWrapper(row_mat), b, x, Base::m_preconditioner, m_iterations, m_error);
246  m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
247  }
248 
249  protected:
250 };
251 
252 } // end namespace Eigen
253 
254 #endif // EIGEN_MINRES_H
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Namespace containing all symbols from the Eigen library.
MINRES()
Definition: MINRES.h:207
~MINRES()
Definition: MINRES.h:223
A minimal residual solver for sparse symmetric problems.
Definition: MINRES.h:137
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
MINRES(const EigenBase< MatrixDerived > &A)
Definition: MINRES.h:220