Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
FlexibleGMRESSolver.h
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2026 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* solver_fgmres.h (C) 2000-2026 */
9/* */
10/* Flexible GMRES method solver. */
11/*---------------------------------------------------------------------------*/
12#ifndef ARCCORE_ALINA_FGMRESSOLVER_H
13#define ARCCORE_ALINA_FGMRESSOLVER_H
14/*---------------------------------------------------------------------------*/
15/*---------------------------------------------------------------------------*/
16/*
17 * This file is based on the work on AMGCL library (version march 2026)
18 * which can be found at https://github.com/ddemidov/amgcl.
19 *
20 * Copyright (c) 2012-2022 Denis Demidov <dennis.demidov@gmail.com>
21 * SPDX-License-Identifier: MIT
22 */
23/*---------------------------------------------------------------------------*/
24/*---------------------------------------------------------------------------*/
25
26#include "arccore/alina/SolverUtils.h"
27#include "arccore/alina/SolverBase.h"
28
29/*---------------------------------------------------------------------------*/
30/*---------------------------------------------------------------------------*/
31
32namespace Arcane::Alina
33{
34
35/*---------------------------------------------------------------------------*/
36/*---------------------------------------------------------------------------*/
40struct FlexibleGMRESSolverParams
41{
42 using params = FlexibleGMRESSolverParams;
43
45 Int32 M = 30;
46
49
51 double tol = 1.0e-8;
52
54 double abstol = std::numeric_limits<double>::min();
55
61 bool ns_search = false;
62
64 bool verbose = false;
65
66 FlexibleGMRESSolverParams() = default;
67
68 FlexibleGMRESSolverParams(const PropertyTree& p)
69 : ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, M)
70 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, maxiter)
71 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, tol)
72 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, abstol)
73 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, ns_search)
74 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, verbose)
75 {
76 p.check_params({ "M", "maxiter", "tol", "abstol", "ns_search", "verbose" });
77 }
78
79 void get(PropertyTree& p, const std::string& path) const
80 {
81 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, M);
82 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, maxiter);
83 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, tol);
84 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, abstol);
85 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, ns_search);
86 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, verbose);
87 }
88};
89
90/*---------------------------------------------------------------------------*/
91/*---------------------------------------------------------------------------*/
98template <class Backend, class InnerProduct = detail::default_inner_product>
100: public SolverBase
101{
102 public:
103
104 typedef Backend backend_type;
105
106 typedef typename Backend::vector vector;
107 typedef typename Backend::value_type value_type;
108 typedef typename Backend::params backend_params;
109
110 typedef typename math::scalar_of<value_type>::type scalar_type;
111 typedef typename math::rhs_of<value_type>::type rhs_type;
112 typedef typename math::inner_product_impl<rhs_type>::return_type coef_type;
113
114 using params = FlexibleGMRESSolverParams;
115
118 const params& prm = params(),
119 const backend_params& bprm = backend_params(),
120 const InnerProduct& inner_product = InnerProduct())
121 : prm(prm)
122 , n(n)
123 , H(prm.M + 1, prm.M)
124 , s(prm.M + 1)
125 , cs(prm.M + 1)
126 , sn(prm.M + 1)
127 , r(Backend::create_vector(n, bprm))
128 , inner_product(inner_product)
129 {
130 v.reserve(prm.M + 1);
131 for (unsigned i = 0; i <= prm.M; ++i)
132 v.push_back(Backend::create_vector(n, bprm));
133
134 z.reserve(prm.M);
135 for (unsigned i = 0; i < prm.M; ++i)
136 z.push_back(Backend::create_vector(n, bprm));
137 }
138
154 template <class Matrix, class Precond, class Vec1, class Vec2>
155 SolverResult operator()(Matrix const& A, Precond const& P, Vec1 const& rhs, Vec2& x) const
156 {
157 ScopedStreamModifier ss(std::cout);
158
159 scalar_type norm_rhs = norm(rhs);
160 if (norm_rhs < Alina::detail::eps<scalar_type>(1)) {
161 if (prm.ns_search) {
162 norm_rhs = math::identity<scalar_type>();
163 }
164 else {
165 backend::clear(x);
166 return SolverResult(0, norm_rhs);
167 }
168 }
169
170 scalar_type eps = std::max(prm.tol * norm_rhs, prm.abstol);
171 scalar_type norm_r = math::zero<scalar_type>();
172
173 unsigned iter = 0;
174 while (true) {
175 backend::residual(rhs, A, x, *v[0]);
176
177 // -- Check stopping condition
178 if ((norm_r = norm(*v[0])) < eps || iter >= prm.maxiter)
179 break;
180
181 // -- Inner GMRES iteration
182 std::fill(s.begin(), s.end(), 0);
183 s[0] = norm_r;
184
185 backend::axpby(math::inverse(norm_r), *v[0], math::zero<scalar_type>(), *v[0]);
186
187 unsigned j = 0;
188 while (true) {
189 // -- Arnoldi process
190 //
191 // Build an orthonormal basis V and matrix H such that
192 // A V_{i-1} = V_{i} H
193
194 vector& v_new = *v[j + 1];
195
196 P.apply(*v[j], *z[j]);
197 backend::spmv(math::identity<scalar_type>(), A, *z[j],
198 math::zero<scalar_type>(), v_new);
199
200 for (unsigned k = 0; k <= j; ++k) {
201 H(k, j) = inner_product(v_new, *v[k]);
202 backend::axpby(-H(k, j), *v[k], math::identity<scalar_type>(), v_new);
203 }
204 H(j + 1, j) = norm(v_new);
205
206 backend::axpby(math::inverse(H(j + 1, j)), v_new, math::zero<scalar_type>(), v_new);
207
208 for (unsigned k = 0; k < j; ++k)
209 detail::apply_plane_rotation(H(k, j), H(k + 1, j), cs[k], sn[k]);
210
211 detail::generate_plane_rotation(H(j, j), H(j + 1, j), cs[j], sn[j]);
212 detail::apply_plane_rotation(H(j, j), H(j + 1, j), cs[j], sn[j]);
213 detail::apply_plane_rotation(s[j], s[j + 1], cs[j], sn[j]);
214
215 scalar_type inner_res = std::abs(s[j + 1]);
216
217 if (prm.verbose && iter % 5 == 0)
218 std::cout << iter << "\t" << std::scientific << inner_res / norm_rhs << std::endl;
219
220 // Check for termination
221 ++j, ++iter;
222 if (iter >= prm.maxiter || j >= prm.M || inner_res <= eps)
223 break;
224 }
225
226 // -- GMRES terminated: eval solution
227 for (unsigned i = j; i-- > 0;) {
228 s[i] /= H(i, i);
229 for (unsigned k = 0; k < i; ++k)
230 s[k] -= H(k, i) * s[i];
231 }
232
233 backend::lin_comb(j, s, z, math::identity<scalar_type>(), x);
234 }
235
236 return SolverResult(iter, norm_r / norm_rhs);
237 }
238
249 template <class Precond, class Vec1, class Vec2>
250 SolverResult operator()(Precond const& P, Vec1 const& rhs, Vec2& x) const
251 {
252 return (*this)(P.system_matrix(), P, rhs, x);
253 }
254
255 size_t bytes() const
256 {
257 size_t b = 0;
258
259 b += H.size() * sizeof(coef_type);
260 b += backend::bytes(s);
261 b += backend::bytes(cs);
262 b += backend::bytes(sn);
263 b += backend::bytes(*r);
264
265 for (const auto& x : v)
266 b += backend::bytes(*x);
267 for (const auto& x : z)
268 b += backend::bytes(*x);
269
270 return b;
271 }
272
273 friend std::ostream& operator<<(std::ostream& os, const FlexibleGMRESSolver& s)
274 {
275 return os << "Type: FGMRES(" << s.prm.M << ")"
276 << "\nUnknowns: " << s.n
277 << "\nMemory footprint: " << human_readable_memory(s.bytes())
278 << "\n";
279 }
280
281 public:
282
283 params prm;
284
285 private:
286
287 size_t n;
288
289 mutable multi_array<coef_type, 2> H;
290 mutable std::vector<coef_type> s, cs, sn;
291 std::shared_ptr<vector> r;
292 std::vector<std::shared_ptr<vector>> v;
293 std::vector<std::shared_ptr<vector>> z;
294
295 InnerProduct inner_product;
296
297 template <class Vec>
298 scalar_type norm(const Vec& x) const
299 {
300 return std::abs(sqrt(inner_product(x, x)));
301 }
302};
303
304/*---------------------------------------------------------------------------*/
305/*---------------------------------------------------------------------------*/
306
307} // namespace Arcane::Alina
308
309/*---------------------------------------------------------------------------*/
310/*---------------------------------------------------------------------------*/
311
312#endif
Flexible GMRES method. \rst Flexible version of the GMRES method [Saad03]_. \endrst.
SolverResult operator()(Matrix const &A, Precond const &P, Vec1 const &rhs, Vec2 &x) const
Computes the solution for the given system matrix.
FlexibleGMRESSolver(size_t n, const params &prm=params(), const backend_params &bprm=backend_params(), const InnerProduct &inner_product=InnerProduct())
Preallocates necessary data structures for the system of size n.
size_t bytes() const
Memory used in bytes.
SolverResult operator()(Precond const &P, Vec1 const &rhs, Vec2 &x) const
Computes the solution for the given right-hand side.
Save ostream flags in constructor, restore in destructor.
Base class for solvers.
Definition SolverBase.h:31
Result of a solving.
Definition AlinaUtils.h:52
Matrix class, to be used by user.
std::int32_t Int32
Type entier signé sur 32 bits.
Parameters for FlexibleGMRES solver.
bool ns_search
Ignore the trivial solution x=0 when rhs is zero.
double tol
Target relative residual error.
bool verbose
Verbose output (show iterations and error)
Int32 maxiter
Maximum number of iterations.
Int32 M
Number of inner GMRES iterations per each outer iteration.
double abstol
Target absolute residual error.