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