GetFEM  5.5
gmm_solver_cg.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 cg.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_cg.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 Conjugate gradient iterative solver.
69 */
70 #ifndef GMM_SOLVER_CG_H__
71 #define GMM_SOLVER_CG_H__
72 
73 #include "gmm_kernel.h"
74 #include "gmm_iter.h"
75 
76 namespace gmm {
77 
78  /* ******************************************************************** */
79  /* conjugate gradient */
80  /* (preconditionned, with parametrable additional scalar product) */
81  /* ******************************************************************** */
82 
83  template <typename Matrix, typename Matps, typename Precond,
84  typename Vector1, typename Vector2>
85  void cg(const Matrix& A, Vector1& x, const Vector2& b, const Matps& PS,
86  const Precond &P, iteration &iter) {
87 
88  typedef typename temporary_dense_vector<Vector1>::vector_type temp_vector;
89  typedef typename linalg_traits<Vector1>::value_type T;
90 
91  T rho, rho_1(0), a;
92  temp_vector p(vect_size(x)), q(vect_size(x)), r(vect_size(x)),
93  z(vect_size(x));
94  iter.set_rhsnorm(gmm::sqrt(gmm::abs(vect_hp(PS, b, b))));
95 
96  if (iter.get_rhsnorm() == 0.0)
97  clear(x);
98  else {
99  mult(A, scaled(x, T(-1)), b, r);
100  mult(P, r, z);
101  rho = vect_hp(PS, z, r);
102  copy(z, p);
103 
104  while (!iter.finished_vect(r)) {
105 
106  if (!iter.first()) {
107  mult(P, r, z);
108  rho = vect_hp(PS, z, r);
109  add(z, scaled(p, rho / rho_1), p);
110  }
111  mult(A, p, q);
112 
113  a = rho / vect_hp(PS, q, p);
114  add(scaled(p, a), x);
115  add(scaled(q, -a), r);
116  rho_1 = rho;
117 
118  ++iter;
119  }
120  }
121  }
122 
123  template <typename Matrix, typename Matps, typename Precond,
124  typename Vector1, typename Vector2>
125  void cg(const Matrix& A, Vector1& x, const Vector2& b, const Matps& PS,
126  const gmm::identity_matrix &, iteration &iter) {
127 
128  typedef typename temporary_dense_vector<Vector1>::vector_type temp_vector;
129  typedef typename linalg_traits<Vector1>::value_type T;
130 
131  T rho, rho_1(0), a;
132  temp_vector p(vect_size(x)), q(vect_size(x)), r(vect_size(x));
133  iter.set_rhsnorm(gmm::sqrt(gmm::abs(vect_hp(PS, b, b))));
134 
135  if (iter.get_rhsnorm() == 0.0)
136  clear(x);
137  else {
138  mult(A, scaled(x, T(-1)), b, r);
139  rho = vect_hp(PS, r, r);
140  copy(r, p);
141 
142  while (!iter.finished_vect(r)) {
143 
144  if (!iter.first()) {
145  rho = vect_hp(PS, r, r);
146  add(r, scaled(p, rho / rho_1), p);
147  }
148  mult(A, p, q);
149  a = rho / vect_hp(PS, q, p);
150  add(scaled(p, a), x);
151  add(scaled(q, -a), r);
152  rho_1 = rho;
153  ++iter;
154  }
155  }
156  }
157 
158  template <typename Matrix, typename Matps, typename Precond,
159  typename Vector1, typename Vector2> inline
160  void cg(const Matrix& A, const Vector1& x, const Vector2& b, const Matps& PS,
161  const Precond &P, iteration &iter)
162  { cg(A, linalg_const_cast(x), b, PS, P, iter); }
163 
164  template <typename Matrix, typename Precond,
165  typename Vector1, typename Vector2> inline
166  void cg(const Matrix& A, Vector1& x, const Vector2& b,
167  const Precond &P, iteration &iter)
168  { cg(A, x , b, identity_matrix(), P, iter); }
169 
170  template <typename Matrix, typename Precond,
171  typename Vector1, typename Vector2> inline
172  void cg(const Matrix& A, const Vector1& x, const Vector2& b,
173  const Precond &P, iteration &iter)
174  { cg(A, x , b , identity_matrix(), P , iter); }
175 
176 }
177 
178 
179 #endif // GMM_SOLVER_CG_H__
void copy(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:976
strongest_value_type< V1, V2 >::value_type vect_hp(const V1 &v1, const V2 &v2)
*‍/
Definition: gmm_blas.h:510
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
void add(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:1275
Iteration object.
Include the base gmm files.