1
// -*- c++ -*- (enables emacs c++ mode)
2
//===========================================================================
4
// Copyright (C) 1997-2008 Yves Renard
6
// This file is a part of GETFEM++
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 2.1 of the License, or
11
// (at your option) any later version.
12
// This program is distributed in the hope that it will be useful, but
13
// WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
14
// or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
15
// License for more details.
16
// You should have received a copy of the GNU Lesser General Public License
17
// along with this program; if not, write to the Free Software Foundation,
18
// Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.
20
// As a special exception, you may use this file as part of a free software
21
// library without restriction. Specifically, if other files instantiate
22
// templates or use macros or inline functions from this file, or you compile
23
// this file and link it with other files to produce an executable, this
24
// file does not by itself cause the resulting executable to be covered by
25
// the GNU General Public License. This exception does not however
26
// invalidate any other reasons why the executable file might be covered by
27
// the GNU General Public License.
29
//===========================================================================
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
//===========================================================================
36
// Copyright (c) 1997-2001, The Trustees of Indiana University.
37
// 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:
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 California, Berkeley 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.
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.
62
//===========================================================================
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.
71
#ifndef GMM_KRYLOV_GMRES_H
72
#define GMM_KRYLOV_GMRES_H
74
#include "gmm_kernel.h"
76
#include "gmm_modified_gram_schmidt.h"
80
/** Generalized Minimum Residual
82
This solve the unsymmetric linear system Ax = b using restarted GMRES.
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
88
template <typename Mat, typename Vec, typename VecB, typename Precond,
90
void gmres(const Mat &A, Vec &x, const VecB &b, const Precond &M,
91
int restart, iteration &outer, Basis& KS) {
93
typedef typename linalg_traits<Vec>::value_type T;
94
typedef typename number_traits<T>::magnitude_type R;
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);
100
double t_ref, t_prec = MPI_Wtime(), t_tot = 0;
101
static double tmult_tot = 0.0;
103
cout << "GMRES " << endl;
106
outer.set_rhsnorm(gmm::vect_norm2(r));
107
if (outer.get_rhsnorm() == 0.0) { clear(x); return; }
109
mult(A, scaled(x, -1), b, w);
111
R beta = gmm::vect_norm2(r), beta_old = beta;
114
iteration inner = outer;
115
inner.reduce_noisy();
116
inner.set_maxiter(restart);
117
inner.set_name("GMRes inner");
119
while (! outer.finished(beta)) {
121
gmm::copy(gmm::scaled(r, R(1)/beta), KS[0]);
125
size_type i = 0; inner.init();
130
orthogonalize(KS, mat_col(H, i), i);
131
R a = gmm::vect_norm2(KS[i+1]);
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]);
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]);
141
++inner, ++outer, ++i;
142
} while (! inner.finished(gmm::abs(s[i])));
144
upper_tri_solve(H, s, i, false);
145
combine(KS, s, x, i);
146
mult(A, gmm::scaled(x, -1), b, w);
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;
152
if (outer.get_noisy()) cout << "Gmres is blocked, exiting\n";
156
t_tot = MPI_Wtime() - t_ref;
157
cout << "temps GMRES : " << t_tot << endl;
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);