GetFEM  5.4.3
gmm_solver_bicgstab.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 bicgstab.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_bicgstab.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 BiCGStab iterative solver.
70 */
71 
72 #ifndef GMM_SOLVER_BICGSTAB_H__
73 #define GMM_SOLVER_BICGSTAB_H__
74 
75 #include "gmm_kernel.h"
76 #include "gmm_iter.h"
77 
78 namespace gmm {
79 
80  /* ******************************************************************** */
81  /* BiConjugate Gradient Stabilized */
82  /* (preconditionned, with parametrable scalar product) */
83  /* ******************************************************************** */
84 
85  template <typename Matrix, typename Vector, typename VectorB,
86  typename Preconditioner>
87  void bicgstab(const Matrix& A, Vector& x, const VectorB& b,
88  const Preconditioner& M, iteration &iter) {
89 
90  typedef typename linalg_traits<Vector>::value_type T;
91  typedef typename number_traits<T>::magnitude_type R;
92  typedef typename temporary_dense_vector<Vector>::vector_type temp_vector;
93 
94  T rho_1, rho_2(0), alpha(0), beta, omega(0);
95  temp_vector p(vect_size(x)), phat(vect_size(x)), s(vect_size(x)),
96  shat(vect_size(x)),
97  t(vect_size(x)), v(vect_size(x)), r(vect_size(x)), rtilde(vect_size(x));
98 
99  gmm::mult(A, gmm::scaled(x, -T(1)), b, r);
100  gmm::copy(r, rtilde);
101  R norm_r = gmm::vect_norm2(r);
102  iter.set_rhsnorm(gmm::vect_norm2(b));
103 
104  if (iter.get_rhsnorm() == 0.0) { clear(x); return; }
105 
106  while (!iter.finished(norm_r)) {
107 
108  rho_1 = gmm::vect_sp(rtilde, r);
109  if (rho_1 == T(0)) {
110  if (iter.get_maxiter() == size_type(-1))
111  { GMM_ASSERT1(false, "Bicgstab failed to converge"); }
112  else { GMM_WARNING1("Bicgstab failed to converge"); return; }
113  }
114 
115  if (iter.first())
116  gmm::copy(r, p);
117  else {
118  if (omega == T(0)) {
119  if (iter.get_maxiter() == size_type(-1))
120  { GMM_ASSERT1(false, "Bicgstab failed to converge"); }
121  else { GMM_WARNING1("Bicgstab failed to converge"); return; }
122  }
123 
124  beta = (rho_1 / rho_2) * (alpha / omega);
125 
126  gmm::add(gmm::scaled(v, -omega), p);
127  gmm::add(r, gmm::scaled(p, beta), p);
128  }
129  gmm::mult(M, p, phat);
130  gmm::mult(A, phat, v);
131  alpha = rho_1 / gmm::vect_sp(v, rtilde);
132  gmm::add(r, gmm::scaled(v, -alpha), s);
133 
134  if (iter.finished_vect(s))
135  { gmm::add(gmm::scaled(phat, alpha), x); break; }
136 
137  gmm::mult(M, s, shat);
138  gmm::mult(A, shat, t);
139  omega = gmm::vect_sp(t, s) / gmm::vect_norm2_sqr(t);
140 
141  gmm::add(gmm::scaled(phat, alpha), x);
142  gmm::add(gmm::scaled(shat, omega), x);
143  gmm::add(s, gmm::scaled(t, -omega), r);
144  norm_r = gmm::vect_norm2(r);
145  rho_2 = rho_1;
146 
147  ++iter;
148  }
149  }
150 
151  template <typename Matrix, typename Vector, typename VectorB,
152  typename Preconditioner>
153  void bicgstab(const Matrix& A, const Vector& x, const VectorB& b,
154  const Preconditioner& M, iteration &iter)
155  { bicgstab(A, linalg_const_cast(x), b, M, iter); }
156 
157 }
158 
159 
160 #endif // GMM_SOLVER_BICGSTAB_H__
void copy(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:978
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
Definition: gmm_blas.h:558
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:545
void clear(L &l)
clear (fill with zeros) a vector or matrix.
Definition: gmm_blas.h:59
void mult(const L1 &l1, const L2 &l2, L3 &l3)
*‍/
Definition: gmm_blas.h:1664
strongest_value_type< V1, V2 >::value_type vect_sp(const V1 &v1, const V2 &v2)
*‍/
Definition: gmm_blas.h:264
void add(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:1277
Iteration object.
Include the base gmm files.
size_t size_type
used as the common size type in the library
Definition: bgeot_poly.h:49
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:47