GetFEM  5.5
gmm_solver_bicgstab.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 bicgstab.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_bicgstab.h
64  @author Andrew Lumsdaine <lums@osl.iu.edu>
65  @author Lie-Quan Lee <llee@osl.iu.edu>
66  @author Yves Renard <Yves.Renard@insa-lyon.fr>
67  @date October 13, 2002.
68  @brief BiCGStab iterative solver.
69 */
70 
71 #ifndef GMM_SOLVER_BICGSTAB_H__
72 #define GMM_SOLVER_BICGSTAB_H__
73 
74 #include "gmm_kernel.h"
75 #include "gmm_iter.h"
76 
77 namespace gmm {
78 
79  /* ******************************************************************** */
80  /* BiConjugate Gradient Stabilized */
81  /* (preconditionned, with parametrable scalar product) */
82  /* ******************************************************************** */
83 
84  template <typename Matrix, typename Vector, typename VectorB,
85  typename Preconditioner>
86  void bicgstab(const Matrix& A, Vector& x, const VectorB& b,
87  const Preconditioner& M, iteration &iter) {
88 
89  typedef typename linalg_traits<Vector>::value_type T;
90  typedef typename number_traits<T>::magnitude_type R;
91  typedef typename temporary_dense_vector<Vector>::vector_type temp_vector;
92 
93  T rho_1, rho_2(0), alpha(0), beta, omega(0);
94  temp_vector p(vect_size(x)), phat(vect_size(x)), s(vect_size(x)),
95  shat(vect_size(x)),
96  t(vect_size(x)), v(vect_size(x)), r(vect_size(x)), rtilde(vect_size(x));
97 
98  gmm::mult(A, gmm::scaled(x, -T(1)), b, r);
99  gmm::copy(r, rtilde);
100  R norm_r = gmm::vect_norm2(r);
101  iter.set_rhsnorm(gmm::vect_norm2(b));
102 
103  if (iter.get_rhsnorm() == 0.0) { clear(x); return; }
104 
105  while (!iter.finished(norm_r)) {
106 
107  rho_1 = gmm::vect_sp(rtilde, r);
108  if (rho_1 == T(0)) {
109  if (iter.get_maxiter() == size_type(-1))
110  { GMM_ASSERT1(false, "Bicgstab failed to converge"); }
111  else { GMM_WARNING1("Bicgstab failed to converge"); return; }
112  }
113 
114  if (iter.first())
115  gmm::copy(r, p);
116  else {
117  if (omega == T(0)) {
118  if (iter.get_maxiter() == size_type(-1))
119  { GMM_ASSERT1(false, "Bicgstab failed to converge"); }
120  else { GMM_WARNING1("Bicgstab failed to converge"); return; }
121  }
122 
123  beta = (rho_1 / rho_2) * (alpha / omega);
124 
125  gmm::add(gmm::scaled(v, -omega), p);
126  gmm::add(r, gmm::scaled(p, beta), p);
127  }
128  gmm::mult(M, p, phat);
129  gmm::mult(A, phat, v);
130  alpha = rho_1 / gmm::vect_sp(v, rtilde);
131  gmm::add(r, gmm::scaled(v, -alpha), s);
132 
133  if (iter.finished_vect(s))
134  { gmm::add(gmm::scaled(phat, alpha), x); break; }
135 
136  gmm::mult(M, s, shat);
137  gmm::mult(A, shat, t);
138  omega = gmm::vect_sp(t, s) / gmm::vect_norm2_sqr(t);
139 
140  gmm::add(gmm::scaled(phat, alpha), x);
141  gmm::add(gmm::scaled(shat, omega), x);
142  gmm::add(s, gmm::scaled(t, -omega), r);
143  norm_r = gmm::vect_norm2(r);
144  rho_2 = rho_1;
145 
146  ++iter;
147  }
148  }
149 
150  template <typename Matrix, typename Vector, typename VectorB,
151  typename Preconditioner>
152  void bicgstab(const Matrix& A, const Vector& x, const VectorB& b,
153  const Preconditioner& M, iteration &iter)
154  { bicgstab(A, linalg_const_cast(x), b, M, iter); }
155 
156 }
157 
158 
159 #endif // GMM_SOLVER_BICGSTAB_H__
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
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2_sqr(const V &v)
squared Euclidean norm of a vector.
Definition: gmm_blas.h:543
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
strongest_value_type< V1, V2 >::value_type vect_sp(const V1 &v1, const V2 &v2)
*‍/
Definition: gmm_blas.h:263
void add(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:1275
Iteration object.
Include the base gmm files.
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:48
size_type alpha(short_type n, short_type d)
Return the value of which is the number of monomials of a polynomial of variables and degree .
Definition: bgeot_poly.cc:46