gmres2.hpp 5.31 KB
Newer Older
1
#pragma once
2
3

#include <algorithm>
4

5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
#include <boost/numeric/mtl/concept/collection.hpp>
#include <boost/numeric/mtl/vector/dense_vector.hpp>
#include <boost/numeric/mtl/matrix/dense2D.hpp>
#include <boost/numeric/mtl/matrix/multi_vector.hpp>
#include <boost/numeric/mtl/operation/givens.hpp>
#include <boost/numeric/mtl/operation/two_norm.hpp>
#include <boost/numeric/mtl/utility/exception.hpp>
#include <boost/numeric/mtl/utility/irange.hpp>

#include "details.hpp"

namespace itl
{
  /// Generalized Minimal Residual method (without restart)
  /** It computes at most kmax_in iterations (or size(x) depending on what is smaller)
      regardless on whether the termination criterion is reached or not.   **/
  template <typename Matrix, typename Vector, typename LeftPreconditioner, typename RightPreconditioner, typename Iteration>
  int gmres2_full(const Matrix& A, Vector& x, const Vector& b,
                  LeftPreconditioner& L, RightPreconditioner& R, Iteration& iter)
  {
    using mtl::irange;
    using std::abs;
    using math::reciprocal;
    using mtl::conj;
    typedef typename mtl::Collection<Vector>::value_type Scalar;
    typedef typename mtl::Collection<Vector>::size_type  Size;

    if (size(b) == 0) throw mtl::logic_error("empty rhs vector");

    const Scalar zero= math::zero(Scalar()), breakdown_tol= 1.e-30, kappa = 10.0; // kappa in [1.25, 100]
    Scalar       rho, hr, bnrm2, temp;
    Size         k, kmax(std::min(size(x), Size(iter.max_iterations() - iter.iterations())));
    Vector       w(b - A *x), r(solve(L,w));
38
39
40
    mtl::mat::multi_vector<Vector>   V(Vector(resource(x), zero), kmax+1);
    mtl::vec::dense_vector<Scalar>   sn(kmax, zero), cs(kmax, zero), s(kmax+1, zero), y(kmax, zero);  // replicated in distributed solvers
    mtl::mat::dense2D<Scalar>        H(kmax+1, kmax);
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118

    bnrm2 = two_norm(b);
    if (bnrm2 == zero)
    {
      set_to_zero(x);
      // b == 0 => solution = 0
      return iter.terminate(bnrm2);
    }

    rho = two_norm(r);				// norm of preconditioned residual
    if (iter.finished(rho))			// initial guess is good enough solution
      return iter;

    V.vector(0) = r * reciprocal(rho);
    H = zero;
    s[0] = rho;

    // GMRES iteration
    for (k= 0; k < kmax && !iter.finished(rho); ++k, ++iter)
    {

      r = A * Vector(solve(R, V.vector(k)));
      w = solve(L, r);
      temp = two_norm(w);

      // orth(V, w, false);
      // modified Gram Schmidt method
      for (Size i= 0; i < k+1; i++)
      {
        H[i][k]= dot(V.vector(i), w);
        w -= H[i][k] * V.vector(i);
      }
      H[k+1][k] = two_norm(w);

      // reorthogonalization, only if "heuristic" condition is fulfilled
      // see [...] for details
      if (H[k+1][k] < temp * reciprocal(kappa))
      {
        for (Size i= 0; i < k+1; i++)
        {
          hr = dot(w, V.vector(i));
          H[i][k] += hr;
          w -= hr * V.vector(i);
        }
        temp = two_norm(w);

        if (temp < H[k+1][k] * reciprocal(kappa))
        {
          set_to_zero(w);
          H[k+1][k] = 0.0;
        }
        else
        {
          H[k+1][k] = temp;
        }
      }

      // watch for breakdown
      if (H[k+1][k] > breakdown_tol)
        V.vector(k+1) = w * reciprocal(H[k+1][k]);
      else
        return iter.fail(2, "GMRES: Singular Hessenbergmatrix - nearly hard breakdown");

      // k Given's rotations
      for(Size i= 0; i < k; i++)
      {
        temp      =  cs[i]*H[i][k] + sn[i]*H[i+1][k];
        H[i+1][k] = -conj(sn[i])*H[i][k] + cs[i]*H[i+1][k];
        H[i][k]   =  temp;
      }

      details::rotmat(H[k][k], H[k+1][k], cs[k], sn[k]);

      s[k+1] = -conj(sn[k])*s[k];
      s[k]   = cs[k]*s[k];
      H[k][k] = cs[k]*H[k][k] + sn[k]*H[k+1][k];
      H[k+1][k] = 0.0;

119
120
      using std::abs;
      rho = abs(s[k+1]);
121
122
123
124
125
126
127
128
129
130
131
132
133
    }

    // reduce k, to get regular matrix
    //     while (k > 0 && std::abs(s[k-1]) <= iter.atol()) k--;

    // iteration is finished -> compute x: solve H*y=g as far as rank of H allows
    irange range(k);
    for (; !range.empty(); --range)
    {
      try
      {
        y[range] = upper_trisolve(H[range][range], s[range]);
      }
134
      catch (mtl::matrix_singular const&)
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
      {
        continue;    // if singular then try with sub-matrix
      }
      break;
    }

    if (range.finish() < k)
      std::cerr << "GMRES orhogonalized with " << k << " vectors but matrix singular, can only use "
                << range.finish() << " vectors!\n";
    if (range.empty())
      return iter.fail(3, "GMRES did not find any direction to correct x");

    r = V.vector(range) * y[range];
    w = solve(R, r);
    x += w;

    r = b - A*x;
    return iter.terminate(r);
  }

  /// Generalized Minimal Residual method with restart
  template <typename Matrix, typename Vector, typename LeftPreconditioner,
            typename RightPreconditioner, typename Iteration>
  int gmres2(const Matrix& A, Vector& x, const Vector& b,
             LeftPreconditioner& L, RightPreconditioner& R,
             Iteration& iter, typename mtl::Collection<Vector>::size_type restart)
  {
    do
    {
      Iteration inner(iter);
      inner.set_max_iterations(std::min(int(iter.iterations()+restart), iter.max_iterations()));
      inner.suppress_resume(true);
      gmres2_full(A, x, b, L, R, inner);
      iter.update_progress(inner);
    }
    while (!iter.finished());

    return iter;
  }

175
} // end namespace itl