GetFEM  5.5
gmm_solver_gmres.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 2002-2026 Yves Renard
5 
6  This file is a part of GetFEM
7 
8  GetFEM is free software; you can redistribute it and/or modify it
9  under the terms of the GNU Lesser General Public License as published
10  by the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version along with the GCC Runtime Library
12  Exception either version 3.1 or (at your option) any later version.
13  This program is distributed in the hope that it will be useful, but
14  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
15  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
16  License and GCC Runtime Library Exception for more details.
17  You should have received a copy of the GNU Lesser General Public License
18  along with this program. If not, see https://www.gnu.org/licenses/.
19 
20  As a special exception, you may use this file as it is a part of a free
21  software library without restriction. Specifically, if other files
22  instantiate templates or use macros or inline functions from this file,
23  or you compile this file and link it with other files to produce an
24  executable, this file does not by itself cause the resulting executable
25  to be covered by the GNU Lesser General Public License. This exception
26  does not however invalidate any other reasons why the executable file
27  might be covered by the GNU Lesser General Public License.
28 
29 ===========================================================================*/
30 
31 // This file is a modified version of gmres.h from ITL.
32 // See http://osl.iu.edu/research/itl/
33 // Following the corresponding Copyright notice.
34 //===========================================================================
35 //
36 // Copyright (c) 1998-2020, University of Notre Dame. All rights reserved.
37 // Redistribution and use in source and binary forms, with or without
38 // modification, are permitted provided that the following conditions are met:
39 //
40 // * Redistributions of source code must retain the above copyright
41 // notice, this list of conditions and the following disclaimer.
42 // * Redistributions in binary form must reproduce the above copyright
43 // notice, this list of conditions and the following disclaimer in the
44 // documentation and/or other materials provided with the distribution.
45 // * Neither the name of the University of Notre Dame nor the
46 // names of its contributors may be used to endorse or promote products
47 // derived from this software without specific prior written permission.
48 //
49 // THIS SOFTWARE IS PROVIDED BY THE TRUSTEES OF INDIANA UNIVERSITY AND
50 // CONTRIBUTORS ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING,
51 // BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
52 // FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE TRUSTEES
53 // OF INDIANA UNIVERSITY AND CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
54 // INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT
55 // NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
56 // DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
57 // THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
58 // (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF
59 // THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
60 //
61 //===========================================================================
62 
63 /**@file gmm_solver_gmres.h
64  @author Andrew Lumsdaine <[email protected]>
65  @author Lie-Quan Lee <[email protected]>
66  @author Yves Renard <[email protected]>
67  @date October 13, 2002.
68  @brief GMRES (Generalized Minimum Residual) iterative solver.
69 */
70 #ifndef GMM_KRYLOV_GMRES_H
71 #define GMM_KRYLOV_GMRES_H
72 
73 #include "gmm_kernel.h"
74 #include "gmm_iter.h"
76 
77 namespace gmm {
78 
79  /** Generalized Minimum Residual
80 
81  This solve the unsymmetric linear system Ax = b using restarted GMRES.
82 
83  See: Y. Saad and M. Schulter. GMRES: A generalized minimum residual
84  algorithm for solving nonsysmmetric linear systems, SIAM
85  J. Sci. Statist. Comp. 7(1986), pp, 856-869
86  */
87  template <typename Mat, typename Vec, typename VecB, typename Precond,
88  typename Basis >
89  void gmres(const Mat &A, Vec &x, const VecB &b, const Precond &M,
90  int restart, iteration &outer, Basis& KS) {
91 
92  typedef typename linalg_traits<Vec>::value_type T;
93  typedef typename number_traits<T>::magnitude_type R;
94 
95  std::vector<T> w(vect_size(x)), r(vect_size(x)), u(vect_size(x));
96  std::vector<T> c_rot(restart+1), s_rot(restart+1), s(restart+1);
97  gmm::dense_matrix<T> H(restart+1, restart);
98 #ifdef GMM_USES_MPI
99  double t_ref, t_prec = MPI_Wtime(), t_tot = 0;
100  static double tmult_tot = 0.0;
101 t_ref = MPI_Wtime();
102  cout << "GMRES " << endl;
103 #endif
104  mult(M,b,r);
105  outer.set_rhsnorm(gmm::vect_norm2(r));
106  if (outer.get_rhsnorm() == 0.0) { clear(x); return; }
107 
108  mult(A, scaled(x, T(-1)), b, w);
109  mult(M, w, r);
110  R beta = gmm::vect_norm2(r), beta_old = beta;
111  int blocked = 0;
112 
113  iteration inner = outer;
114  inner.reduce_noisy();
115  inner.set_maxiter(restart);
116  inner.set_name("GMRes inner");
117 
118  while (! outer.finished(beta)) {
119 
120  gmm::copy(gmm::scaled(r, R(1)/beta), KS[0]);
121  gmm::clear(s);
122  s[0] = beta;
123 
124  size_type i = 0; inner.init();
125 
126  do {
127  mult(A, KS[i], u);
128  mult(M, u, KS[i+1]);
129  orthogonalize(KS, mat_col(H, i), i);
130  R a = gmm::vect_norm2(KS[i+1]);
131  H(i+1, i) = T(a);
132  gmm::scale(KS[i+1], T(1) / a);
133  for (size_type k = 0; k < i; ++k)
134  Apply_Givens_rotation_left(H(k,i), H(k+1,i), c_rot[k], s_rot[k]);
135 
136  Givens_rotation(H(i,i), H(i+1,i), c_rot[i], s_rot[i]);
137  Apply_Givens_rotation_left(H(i,i), H(i+1,i), c_rot[i], s_rot[i]);
138  Apply_Givens_rotation_left(s[i], s[i+1], c_rot[i], s_rot[i]);
139 
140  ++inner, ++outer, ++i;
141  } while (! inner.finished(gmm::abs(s[i])));
142 
143  upper_tri_solve(H, s, i, false);
144  combine(KS, s, x, i);
145  mult(A, gmm::scaled(x, T(-1)), b, w);
146  mult(M, w, r);
147  beta_old = std::min(beta, beta_old); beta = gmm::vect_norm2(r);
148  if (int(inner.get_iteration()) < restart -1 || beta_old <= beta)
149  ++blocked; else blocked = 0;
150  if (blocked > 10) {
151  if (outer.get_noisy()) cout << "Gmres is blocked, exiting\n";
152  break;
153  }
154 #ifdef GMM_USES_MPI
155  t_tot = MPI_Wtime() - t_ref;
156  cout << "temps GMRES : " << t_tot << endl;
157 #endif
158  }
159  }
160 
161 
162  template <typename Mat, typename Vec, typename VecB, typename Precond >
163  void gmres(const Mat &A, Vec &x, const VecB &b,
164  const Precond &M, int restart, iteration& outer) {
165  typedef typename linalg_traits<Vec>::value_type T;
166  modified_gram_schmidt<T> orth(restart, vect_size(x));
167  gmres(A, x, b, M, restart, outer, orth);
168  }
169 
170 }
171 
172 #endif
The Iteration object calculates whether the solution has reached the desired accuracy,...
Definition: gmm_iter.h:52
void copy(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:976
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
Definition: gmm_blas.h:556
void clear(L &l)
clear (fill with zeros) a vector or matrix.
Definition: gmm_blas.h:58
void mult(const L1 &l1, const L2 &l2, L3 &l3)
*‍/
Definition: gmm_blas.h:1663
Iteration object.
Include the base gmm files.
Modified Gram-Schmidt orthogonalization.
void gmres(const Mat &A, Vec &x, const VecB &b, const Precond &M, int restart, iteration &outer, Basis &KS)
Generalized Minimum Residual.