Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
BiCGStabLSolver.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_bicgstabl.h (C) 2000-2026 */
9/* */
10/* BiCGStab(L) iterative method. . */
11/*---------------------------------------------------------------------------*/
12#ifndef ARCCORE_ALINA_BICGSTABL_H
13#define ARCCORE_ALINA_BICGSTABL_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 *
27 * The code is ported from PETSC BCGSL [1] and is based on [2].
28 *
29 * [1] http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPBCGSL.html
30 * [2] Fokkema, Diederik R. Enhanced implementation of BiCGstab (l) for solving
31 * linear systems of equations. Universiteit Utrecht. Mathematisch Instituut,
32 * 1996.
33
34 The original code came with the following license:
35
36 Copyright (c) 1991-2014, UChicago Argonne, LLC and the PETSc Development Team
37 All rights reserved.
38
39 Redistribution and use in source and binary forms, with or without modification,
40 are permitted provided that the following conditions are met:
41
42 * Redistributions of source code must retain the above copyright notice, this
43 list of conditions and the following disclaimer.
44 * Redistributions in binary form must reproduce the above copyright notice, this
45 list of conditions and the following disclaimer in the documentation and/or
46 other materials provided with the distribution.
47
48 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
49 ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
50 WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
51 DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR
52 ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
53 (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
54 LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
55 ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
56 (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
57 SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
58*/
59/*---------------------------------------------------------------------------*/
60/*---------------------------------------------------------------------------*/
61
62#include "arccore/alina/ValueTypeInterface.h"
63#include "arccore/alina/SolverUtils.h"
64#include "arccore/alina/QRFactorizationImpl.h"
65#include "arccore/alina/SolverBase.h"
66
67/*---------------------------------------------------------------------------*/
68/*---------------------------------------------------------------------------*/
69
70namespace Arcane::Alina
71{
72
73/*---------------------------------------------------------------------------*/
74/*---------------------------------------------------------------------------*/
78struct BiCGStabLSolverParams
79{
80 using params = BiCGStabLSolverParams;
81
82 // Order of the method.
83 int L = 2;
84
85 // Threshold used to decide when to refresh computed residuals.
86 double delta = 0.0;
87
88 // Use a convex function of the MinRes and OR polynomials
89 // after the BiCG step instead of default MinRes
90 bool convex = true;
91
92 // Preconditioning kind (left/right).
93 ePreconditionerSideType pside = ePreconditionerSideType::right;
94
95 // Maximum number of iterations.
96 Int32 maxiter = 100;
97
98 // Target relative residual error.
99 double tol = 1e-8;
100
101 // Target absolute residual error.
102 double abstol = std::numeric_limits<double>::min();
103
109 bool ns_search = false;
110
112 bool verbose = false;
113
114 BiCGStabLSolverParams() = default;
115
116 BiCGStabLSolverParams(const PropertyTree& p)
117 : ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, L)
118 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, delta)
119 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, convex)
120 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, pside)
121 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, maxiter)
122 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, tol)
123 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, abstol)
124 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, ns_search)
125 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, verbose)
126 {
127 p.check_params({ "L", "delta", "convex", "pside", "maxiter", "tol", "abstol", "ns_search", "verbose" });
128 }
129
130 void get(PropertyTree& p, const std::string& path) const
131 {
132 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, L);
133 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, delta);
134 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, convex);
135 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, pside);
136 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, maxiter);
137 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, tol);
138 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, abstol);
139 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, ns_search);
140 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, verbose);
141 }
142};
143
144/*---------------------------------------------------------------------------*/
145/*---------------------------------------------------------------------------*/
149template <class Backend, class InnerProduct = detail::default_inner_product>
151: public SolverBase
152{
153 public:
154
155 using backend_type = Backend;
156 using BackendType = Backend;
157
158 typedef typename Backend::vector vector;
159 typedef typename Backend::value_type value_type;
160 typedef typename Backend::params backend_params;
161
162 typedef typename math::scalar_of<value_type>::type scalar_type;
163
164 typedef typename math::inner_product_impl<
165 typename math::rhs_of<value_type>::type>::return_type coef_type;
166
167 using params = BiCGStabLSolverParams;
168
171 const params& prm = params(),
172 const backend_params& backend_prm = backend_params(),
173 const InnerProduct& inner_product = InnerProduct())
174 : prm(prm)
175 , n(n)
176 , Rt(Backend::create_vector(n, backend_prm))
177 , X(Backend::create_vector(n, backend_prm))
178 , B(Backend::create_vector(n, backend_prm))
179 , T(Backend::create_vector(n, backend_prm))
180 , R(prm.L + 1)
181 , U(prm.L + 1)
182 , MZa(prm.L + 1, prm.L + 1)
183 , MZb(prm.L + 1, prm.L + 1)
184 , Y0(prm.L + 1)
185 , YL(prm.L + 1)
186 , inner_product(inner_product)
187 {
188 precondition(prm.L > 0, "L in BiCGStab(L) should be >=1");
189
190 for (int i = 0; i <= prm.L; ++i) {
191 R[i] = Backend::create_vector(n, backend_prm);
192 U[i] = Backend::create_vector(n, backend_prm);
193 }
194 }
195
196 /*
197 * \brief Computes the solution for the given system matrix.
198 *
199 * Computes the solution for the given system matrix \p A and the
200 * right-hand side \p rhs. Returns the number of iterations made and
201 * the achieved residual as a ``std::tuple``. The solution vector
202 * \p x provides initial approximation in input and holds the computed
203 * solution on output.
204 *
205 * The system matrix may differ from the matrix used during
206 * initialization. This may be used for the solution of non-stationary
207 * problems with slowly changing coefficients. There is a strong chance
208 * that a preconditioner built for a time step will act as a reasonably
209 * good preconditioner for several subsequent time steps [DeSh12]_.
210 */
211 template <class Matrix, class Precond, class Vec1, class Vec2>
212 SolverResult operator()(const Matrix& A, const Precond& P, const Vec1& rhs, Vec2&& x) const
213 {
214 static const coef_type one = math::identity<coef_type>();
215 static const coef_type zero = math::zero<coef_type>();
216
217 const int L = prm.L;
218
219 ScopedStreamModifier ss(std::cout);
220
221 scalar_type norm_rhs = norm(rhs);
222
223 // Check if there is a trivial solution
224 if (norm_rhs < Alina::detail::eps<scalar_type>(1)) {
225 if (prm.ns_search) {
226 norm_rhs = math::identity<scalar_type>();
227 }
228 else {
229 backend::clear(x);
230 return SolverResult(0, norm_rhs);
231 }
232 }
233
234 if (prm.pside == ePreconditionerSideType::left) {
235 backend::residual(rhs, A, x, *T);
236 P.apply(*T, *B);
237 }
238 else {
239 backend::residual(rhs, A, x, *B);
240 }
241
242 scalar_type zeta0 = norm(*B);
243 scalar_type eps = std::max(prm.tol * norm_rhs, prm.abstol);
244
245 coef_type alpha = zero;
246 coef_type rho0 = one;
247 coef_type omega = one;
248
249 // Go
250 backend::copy(*B, *R[0]);
251 backend::copy(*B, *Rt);
252 backend::clear(*X);
253 backend::clear(*U[0]);
254
255 scalar_type zeta = zeta0;
256 scalar_type rnmax_computed = zeta0;
257 scalar_type rnmax_true = zeta0;
258
259 size_t iter = 0;
260 for (; iter < prm.maxiter && zeta >= eps; iter += L) {
261 // BiCG part
262 rho0 = -omega * rho0;
263
264 for (int j = 0; j < L; ++j) {
265 coef_type rho1 = inner_product(*R[j], *Rt);
266 precondition(!math::is_zero(rho1), "BiCGStab(L) breakdown: diverged (zero rho)");
267
268 coef_type beta = alpha * (rho1 / rho0);
269 rho0 = rho1;
270
271 for (int i = 0; i <= j; ++i)
272 backend::axpby(one, *R[i], -beta, *U[i]);
273
274 preconditioner_spmv(prm.pside, P, A, *U[j], *U[j + 1], *T);
275
276 coef_type sigma = inner_product(*U[j + 1], *Rt);
277 precondition(!math::is_zero(sigma), "BiCGStab(L) breakdown: diverged (zero sigma)");
278 alpha = rho1 / sigma;
279
280 backend::axpby(alpha, *U[0], one, *X);
281
282 for (int i = 0; i <= j; ++i)
283 backend::axpby(-alpha, *U[i + 1], one, *R[i]);
284
285 preconditioner_spmv(prm.pside, P, A, *R[j], *R[j + 1], *T);
286
287 zeta = norm(*R[0]);
288
289 rnmax_computed = std::max(zeta, rnmax_computed);
290 rnmax_true = std::max(zeta, rnmax_true);
291
292 // Check for early exit
293 if (zeta < eps) {
294 iter += j + 1;
295 goto done;
296 }
297 }
298
299 // Polynomial part
300 for (int i = 0; i <= L; ++i) {
301 for (int j = 0; j <= i; ++j) {
302 MZa(i, j) = inner_product(*R[i], *R[j]);
303 }
304 }
305
306 // Symmetrize MZa
307 for (int i = 0; i <= L; ++i) {
308 for (int j = i + 1; j <= L; ++j) {
309 MZa(i, j) = MZa(j, i) = math::adjoint(MZa(j, i));
310 }
311 }
312
313 std::copy(MZa.data(), MZa.data() + MZa.size(), MZb.data());
314
315 if (prm.convex || L == 1) {
316 Y0[0] = -one;
317
318 qr.solve(L, L, MZa.stride(0), MZa.stride(1),
319 &MZa(1, 1), &MZb(0, 1), &Y0[1]);
320 }
321 else {
322 Y0[0] = -one;
323 Y0[L] = zero;
324 qr.solve(L - 1, L - 1, MZa.stride(0), MZa.stride(1),
325 &MZa(1, 1), &MZb(0, 1), &Y0[1]);
326
327 YL[0] = zero;
328 YL[L] = -one;
329 qr.solve(L - 1, L - 1, MZa.stride(0), MZa.stride(1),
330 &MZa(1, 1), &MZb(L, 1), &YL[1], /*computed=*/true);
331
332 coef_type dot0 = zero;
333 coef_type dot1 = zero;
334 coef_type dotA = zero;
335 for (int i = 0; i <= L; ++i) {
336 coef_type s0 = zero;
337 coef_type sL = zero;
338
339 for (int j = 0; j <= L; ++j) {
340 coef_type M = MZb(i, j);
341 s0 += M * Y0[j];
342 sL += M * YL[j];
343 }
344
345 dot0 += Y0[i] * s0;
346 dotA += YL[i] * s0;
347 dot1 += YL[i] * sL;
348 }
349
350 scalar_type kappa0 = sqrt(std::abs(std::real(dot0)));
351 scalar_type kappa1 = sqrt(std::abs(std::real(dot1)));
352 scalar_type kappaA = std::real(dotA);
353
354 if (!math::is_zero(kappa0) && !math::is_zero(kappa1)) {
355 scalar_type ghat;
356 if (kappaA < 0.7 * kappa0 * kappa1) {
357 ghat = (kappaA < 0) ? -0.7 * kappa0 / kappa1 : 0.7 * kappa0 / kappa1;
358 }
359 else {
360 ghat = kappaA / (kappa1 * kappa1);
361 }
362
363 for (int i = 0; i <= L; ++i)
364 Y0[i] -= ghat * YL[i];
365 }
366 }
367
368 omega = Y0[L];
369 for (int h = L; h > 0 && math::is_zero(omega); --h)
370 omega = Y0[h];
371 precondition(!math::is_zero(omega), "BiCGStab(L) breakdown: diverged (zero omega)");
372
373 backend::lin_comb(L, &Y0[1], &R[0], one, *X);
374
375 for (int i = 1; i <= L; ++i)
376 Y0[i] = -one * Y0[i];
377
378 backend::lin_comb(L, &Y0[1], &U[1], one, *U[0]);
379 backend::lin_comb(L, &Y0[1], &R[1], one, *R[0]);
380
381 for (int i = 1; i <= L; ++i)
382 Y0[i] = -one * Y0[i];
383
384 zeta = norm(*R[0]);
385
386 // Accurate update
387 if (prm.delta > 0) {
388 rnmax_computed = std::max(zeta, rnmax_computed);
389 rnmax_true = std::max(zeta, rnmax_true);
390
391 bool update_x = zeta < prm.delta * zeta0 && zeta0 <= rnmax_computed;
392
393 if ((zeta < prm.delta * rnmax_true && zeta <= rnmax_true) || update_x) {
394 preconditioner_spmv(prm.pside, P, A, *X, *R[0], *T);
395 backend::axpby(one, *B, -one, *R[0]);
396 rnmax_true = zeta;
397
398 if (update_x) {
399 if (prm.pside == ePreconditionerSideType::left) {
400 backend::axpby(one, *X, one, x);
401 }
402 else {
403 backend::axpby(one, *T, one, x);
404 }
405 backend::clear(*X);
406 backend::copy(*R[0], *B);
407
408 rnmax_computed = zeta;
409 }
410 }
411 }
412 if (prm.verbose && iter % 5 == 0)
413 std::cout << iter << "\t" << std::scientific << zeta / norm_rhs << std::endl;
414 }
415
416 done:
417 if (prm.pside == ePreconditionerSideType::left) {
418 backend::axpby(one, *X, one, x);
419 }
420 else {
421 P.apply(*X, *T);
422 backend::axpby(one, *T, one, x);
423 }
424
425 return SolverResult(iter, zeta / norm_rhs);
426 }
427
438 template <class Precond, class Vec1, class Vec2>
439 SolverResult operator()(const Precond& P, const Vec1& rhs, Vec2&& x) const
440 {
441 return (*this)(P.system_matrix(), P, rhs, x);
442 }
443
444 size_t bytes() const
445 {
446 size_t b = 0;
447
448 b += backend::bytes(*Rt);
449 b += backend::bytes(*X);
450 b += backend::bytes(*B);
451 b += backend::bytes(*T);
452
453 for (const auto& v : R)
454 b += backend::bytes(*v);
455 for (const auto& v : U)
456 b += backend::bytes(*v);
457
458 b += MZa.size() * sizeof(coef_type);
459 b += MZb.size() * sizeof(coef_type);
460
461 b += backend::bytes(Y0);
462 b += backend::bytes(YL);
463
464 b += qr.bytes();
465
466 return b;
467 }
468
469 friend std::ostream& operator<<(std::ostream& os, const BiCGStabLSolver& s)
470 {
471 return os << "Type: BiCGStab(" << s.prm.L << ")"
472 << "\nUnknowns: " << s.n
473 << "\nMemory footprint: " << human_readable_memory(s.bytes())
474 << std::endl;
475 }
476
477 public:
478
479 params prm;
480
481 private:
482
483 size_t n;
484
485 mutable std::shared_ptr<vector> Rt;
486 mutable std::shared_ptr<vector> X;
487 mutable std::shared_ptr<vector> B;
488 mutable std::shared_ptr<vector> T;
489
490 mutable std::vector<std::shared_ptr<vector>> R;
491 mutable std::vector<std::shared_ptr<vector>> U;
492
493 mutable multi_array<coef_type, 2> MZa, MZb;
494 mutable std::vector<coef_type> Y0, YL;
496
497 InnerProduct inner_product;
498
499 template <class Vec>
500 scalar_type norm(const Vec& x) const
501 {
502 return sqrt(math::norm(inner_product(x, x)));
503 }
504};
505
506/*---------------------------------------------------------------------------*/
507/*---------------------------------------------------------------------------*/
508
509} // namespace Arcane::Alina
510
511/*---------------------------------------------------------------------------*/
512/*---------------------------------------------------------------------------*/
513
514#endif
size_t bytes() const
Memory used in bytes.
BiCGStabLSolver(size_t n, const params &prm=params(), const backend_params &backend_prm=backend_params(), const InnerProduct &inner_product=InnerProduct())
Preallocates necessary data structures for the system of size n.
SolverResult operator()(const Precond &P, const Vec1 &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.
__host__ __device__ double sqrt(double v)
Racine carrée de v.
Definition Math.h:135
std::int32_t Int32
Type entier signé sur 32 bits.
Parameters for BiCGStab(L) solver.
bool ns_search
Ignore the trivial solution x=0 when rhs is zero.
bool verbose
Verbose output (show iterations and error)
Default implementation for inner product.