Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
AMG.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/* AMG.h (C) 2000-2026 */
9/* */
10/* An AMG preconditioner. . */
11/*---------------------------------------------------------------------------*/
12#ifndef ARCCORE_ALINA_AMG_H
13#define ARCCORE_ALINA_AMG_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// Temporarily add this pragma because there is a lot of conversion warnings.
27#pragma GCC diagnostic ignored "-Wconversion"
28
29/*---------------------------------------------------------------------------*/
30/*---------------------------------------------------------------------------*/
31
32#include "arccore/alina/BuiltinBackend.h"
33#include "arccore/alina/SolverUtils.h"
34#include "arccore/alina/AlinaUtils.h"
35
36#include <iostream>
37#include <iomanip>
38#include <list>
39#include <memory>
40
41/*---------------------------------------------------------------------------*/
42/*---------------------------------------------------------------------------*/
43
44namespace Arcane::Alina
45{
46
47/*---------------------------------------------------------------------------*/
48/*---------------------------------------------------------------------------*/
67template <class Backend,
68 template <class> class Coarsening,
69 template <class> class Relax>
70class AMG
71{
72 public:
73
74 typedef Backend backend_type;
75
76 typedef typename Backend::value_type value_type;
77 typedef typename Backend::col_type col_type;
78 typedef typename Backend::ptr_type ptr_type;
79 typedef typename Backend::matrix matrix;
80 typedef typename Backend::vector vector;
81
82 typedef Coarsening<Backend> coarsening_type;
83 typedef Relax<Backend> relax_type;
84
85 typedef typename BuiltinBackend<value_type, col_type, ptr_type>::matrix build_matrix;
86
87 typedef typename math::scalar_of<value_type>::type scalar_type;
88
90 typedef typename Backend::params backend_params;
91
98 struct params
99 {
100 typedef coarsening_type::params coarsening_params;
101 typedef relax_type::params relax_params;
102
103 coarsening_params coarsening;
104 relax_params relax;
105
113 Int32 coarse_enough = Backend::direct_solver::coarse_enough();
114
121 bool direct_coarse = true;
122
130 Int32 max_levels = std::numeric_limits<Int32>::max();
131
134
137
140
143
145 bool allow_rebuild = std::is_same<matrix, build_matrix>::value;
146
147 params() = default;
148
149 params(const PropertyTree& p)
150 : ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p, coarsening)
151 , ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p, relax)
152 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, coarse_enough)
153 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, direct_coarse)
154 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, max_levels)
155 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, npre)
156 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, npost)
157 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, ncycle)
158 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, pre_cycles)
159 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, allow_rebuild)
160 {
161 p.check_params( { "coarsening", "relax", "coarse_enough", "direct_coarse", "max_levels", "npre", "npost", "ncycle", "pre_cycles", "allow_rebuild" });
162
163 precondition(max_levels > 0, "max_levels should be positive");
164 }
165
166 void get(PropertyTree& p, const std::string& path = "") const
167 {
168 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, coarsening);
169 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, relax);
170 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, coarse_enough);
171 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, direct_coarse);
172 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, max_levels);
173 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, npre);
174 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, npost);
175 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, ncycle);
176 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, pre_cycles);
177 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, allow_rebuild);
178 }
179 };
180
189 template <class Matrix>
190 explicit AMG(const Matrix& M, const params& p = params(),
191 const backend_params& bprm = backend_params())
192 : prm(p)
193 {
194 auto A = std::make_shared<build_matrix>(M);
195 sort_rows(*A);
196
197 _initialize(A, bprm);
198 }
199
212 explicit AMG(std::shared_ptr<build_matrix> A,
213 const params& p = params(),
214 const backend_params& bprm = backend_params())
215 : prm(p)
216 {
217 _initialize(A, bprm);
218 }
219
226 template <class Matrix>
227 void rebuild(const Matrix& M,
228 const backend_params& bprm = backend_params())
229 {
230 auto A = std::make_shared<build_matrix>(M);
231 sort_rows(*A);
232 rebuild(A, bprm);
233 }
234
241 void rebuild(std::shared_ptr<build_matrix> A,
242 const backend_params& bprm = backend_params())
243 {
244 precondition(prm.allow_rebuild, "allow_rebuild is not set!");
245 precondition(backend::nbRow(*A) == backend::nbRow(system_matrix()) &&
246 backend::nbColumn(*A) == backend::nbRow(*A),
247 "Matrix dimensions differ from the original ones!");
248
249 ARCCORE_ALINA_TIC("rebuild");
250 coarsening_type C(prm.coarsening);
251 for (auto& level : levels) {
252 A = level.rebuild(A, C, prm, bprm);
253 }
254 ARCCORE_ALINA_TOC("rebuild");
255 }
256
263 template <class Vec1, class Vec2>
264 void cycle(const Vec1& rhs, Vec2&& x) const
265 {
266 cycle(levels.begin(), rhs, x);
267 }
268
277 template <class Vec1, class Vec2>
278 void apply(const Vec1& rhs, Vec2&& x) const
279 {
280 if (prm.pre_cycles) {
281 backend::clear(x);
282 for (unsigned i = 0; i < prm.pre_cycles; ++i)
283 cycle(rhs, x);
284 }
285 else {
286 backend::copy(rhs, x);
287 }
288 }
289
291 std::shared_ptr<matrix> system_matrix_ptr() const
292 {
293 return levels.front().A;
294 }
295
296 const matrix& system_matrix() const
297 {
298 return *system_matrix_ptr();
299 }
300
301 size_t bytes() const
302 {
303 size_t b = 0;
304 for (const auto& lvl : levels)
305 b += lvl.bytes();
306 return b;
307 }
308
309 public:
310
311 params prm;
312
313 private:
314
315 struct AMGLevel
316 {
317 size_t m_rows = 0;
318 size_t m_nonzeros = 0;
319
320 std::shared_ptr<vector> f;
321 std::shared_ptr<vector> u;
322 std::shared_ptr<vector> t;
323
324 std::shared_ptr<matrix> A;
325 std::shared_ptr<matrix> P;
326 std::shared_ptr<matrix> R;
327
328 std::shared_ptr<build_matrix> bP;
329 std::shared_ptr<build_matrix> bR;
330
331 std::shared_ptr<typename Backend::direct_solver> solve;
332
333 std::shared_ptr<relax_type> relax;
334
335 size_t bytes() const
336 {
337 size_t b = 0;
338
339 if (f)
340 b += backend::bytes(*f);
341 if (u)
342 b += backend::bytes(*u);
343 if (t)
344 b += backend::bytes(*t);
345
346 if (A)
347 b += backend::bytes(*A);
348 if (P)
349 b += backend::bytes(*P);
350 if (R)
351 b += backend::bytes(*R);
352
353 if (solve)
354 b += backend::bytes(*solve);
355 if (relax)
356 b += backend::bytes(*relax);
357
358 return b;
359 }
360
361 AMGLevel() = default;
362
363 AMGLevel(std::shared_ptr<build_matrix> A, params& prm, const backend_params& bprm)
364 : m_rows(backend::nbRow(*A))
365 , m_nonzeros(backend::nonzeros(*A))
366 {
367 ARCCORE_ALINA_TIC("move to backend");
368 f = Backend::create_vector(m_rows, bprm);
369 u = Backend::create_vector(m_rows, bprm);
370 t = Backend::create_vector(m_rows, bprm);
371 this->A = Backend::copy_matrix(A, bprm);
372 ARCCORE_ALINA_TOC("move to backend");
373
374 ARCCORE_ALINA_TIC("relaxation");
375 relax = std::make_shared<relax_type>(*A, prm.relax, bprm);
376 ARCCORE_ALINA_TOC("relaxation");
377 }
378
379 std::shared_ptr<build_matrix>
380 step_down(std::shared_ptr<build_matrix> A, coarsening_type& C,
381 const backend_params& bprm, bool allow_rebuild)
382 {
383 ARCCORE_ALINA_TIC("transfer operators");
384 std::shared_ptr<build_matrix> P, R;
385
386 try {
387 std::tie(P, R) = C.transfer_operators(*A);
388 }
389 catch (error::empty_level) {
390 ARCCORE_ALINA_TOC("transfer operators");
391 return std::shared_ptr<build_matrix>();
392 }
393
394 sort_rows(*P);
395 sort_rows(*R);
396
397 if (allow_rebuild) {
398 bP = P;
399 bR = R;
400 }
401 ARCCORE_ALINA_TOC("transfer operators");
402
403 ARCCORE_ALINA_TIC("move to backend");
404 this->P = Backend::copy_matrix(P, bprm);
405 this->R = Backend::copy_matrix(R, bprm);
406 ARCCORE_ALINA_TOC("move to backend");
407
408 ARCCORE_ALINA_TIC("coarse operator");
409 A = C.coarse_operator(*A, *P, *R);
410 sort_rows(*A);
411 ARCCORE_ALINA_TOC("coarse operator");
412
413 return A;
414 }
415
416 void create_coarse(std::shared_ptr<build_matrix> A,
417 const backend_params& bprm, bool single_level)
418 {
419 m_rows = backend::nbRow(*A);
420 m_nonzeros = backend::nonzeros(*A);
421
422 u = Backend::create_vector(m_rows, bprm);
423 f = Backend::create_vector(m_rows, bprm);
424
425 solve = Backend::create_solver(A, bprm);
426 if (single_level)
427 this->A = Backend::copy_matrix(A, bprm);
428 }
429
430 std::shared_ptr<build_matrix>
431 rebuild(std::shared_ptr<build_matrix> A,
432 const coarsening_type& C,
433 const params& prm,
434 const backend_params& bprm)
435 {
436 if (this->A) {
437 ARCCORE_ALINA_TIC("move to backend");
438 this->A = Backend::copy_matrix(A, bprm);
439 ARCCORE_ALINA_TOC("move to backend");
440 }
441
442 if (relax) {
443 ARCCORE_ALINA_TIC("relaxation");
444 relax = std::make_shared<relax_type>(*A, prm.relax, bprm);
445 ARCCORE_ALINA_TOC("relaxation");
446 }
447
448 if (solve) {
449 ARCCORE_ALINA_TIC("coarsest level");
450 solve = Backend::create_solver(A, bprm);
451 ARCCORE_ALINA_TOC("coarsest level");
452 }
453
454 if (bP && bR) {
455 ARCCORE_ALINA_TIC("coarse operator");
456 A = C.coarse_operator(*A, *bP, *bR);
457 sort_rows(*A);
458 ARCCORE_ALINA_TOC("coarse operator");
459 }
460
461 return A;
462 }
463
464 size_t rows() const
465 {
466 return m_rows;
467 }
468
469 size_t nonzeros() const
470 {
471 return m_nonzeros;
472 }
473 };
474
475 typedef std::list<AMGLevel>::const_iterator level_iterator;
476
477 std::list<AMGLevel> levels;
478
479 void _initialize(std::shared_ptr<build_matrix> A,
480 const backend_params& bprm = backend_params())
481 {
482 precondition(backend::nbRow(*A) == backend::nbColumn(*A), "Matrix should be square!");
483
484 bool direct_coarse_solve = true;
485
486 coarsening_type C(prm.coarsening);
487 std::cout << "DoInit AMG coarse_enough=" << prm.coarse_enough << "\n";
488 while (backend::nbRow(*A) > prm.coarse_enough) {
489 std::cout << "DoIteration nb_row=" << backend::nbRow(*A) << "\n";
490 levels.push_back(AMGLevel(A, prm, bprm));
491
492 if (levels.size() >= prm.max_levels)
493 break;
494
495 A = levels.back().step_down(A, C, bprm, prm.allow_rebuild);
496 if (!A) {
497 // Zero-sized coarse level. Probably the system matrix on
498 // this level is diagonal, should be easily solvable with a
499 // couple of smoother iterations.
500 direct_coarse_solve = false;
501 break;
502 }
503 }
504
505 if (!A || backend::nbRow(*A) > prm.coarse_enough) {
506 // The coarse matrix is still too big to be solved directly.
507 direct_coarse_solve = false;
508 }
509
510 if (direct_coarse_solve) {
511 ARCCORE_ALINA_TIC("coarsest level");
512 if (prm.direct_coarse) {
513 AMGLevel l;
514 l.create_coarse(A, bprm, levels.empty());
515 levels.push_back(l);
516 }
517 else {
518 levels.push_back(AMGLevel(A, prm, bprm));
519 }
520 ARCCORE_ALINA_TOC("coarsest level");
521 }
522 }
523
524 template <class Vec1, class Vec2>
525 void cycle(level_iterator lvl, const Vec1& rhs, Vec2& x) const
526 {
527 level_iterator nxt = lvl, end = levels.end();
528 ++nxt;
529
530 if (nxt == end) {
531 if (lvl->solve) {
532 ARCCORE_ALINA_TIC("coarse");
533 (*lvl->solve)(rhs, x);
534 ARCCORE_ALINA_TOC("coarse");
535 }
536 else {
537 ARCCORE_ALINA_TIC("relax");
538 for (size_t i = 0; i < prm.npre; ++i)
539 lvl->relax->apply_pre(*lvl->A, rhs, x, *lvl->t);
540 for (size_t i = 0; i < prm.npost; ++i)
541 lvl->relax->apply_post(*lvl->A, rhs, x, *lvl->t);
542 ARCCORE_ALINA_TOC("relax");
543 }
544 }
545 else {
546 for (size_t j = 0; j < prm.ncycle; ++j) {
547 ARCCORE_ALINA_TIC("relax");
548 for (size_t i = 0; i < prm.npre; ++i)
549 lvl->relax->apply_pre(*lvl->A, rhs, x, *lvl->t);
550 ARCCORE_ALINA_TOC("relax");
551
552 backend::residual(rhs, *lvl->A, x, *lvl->t);
553
554 backend::spmv(math::identity<scalar_type>(), *lvl->R, *lvl->t, math::zero<scalar_type>(), *nxt->f);
555
556 backend::clear(*nxt->u);
557 cycle(nxt, *nxt->f, *nxt->u);
558
559 backend::spmv(math::identity<scalar_type>(), *lvl->P, *nxt->u, math::identity<scalar_type>(), x);
560
561 ARCCORE_ALINA_TIC("relax");
562 for (size_t i = 0; i < prm.npost; ++i)
563 lvl->relax->apply_post(*lvl->A, rhs, x, *lvl->t);
564 ARCCORE_ALINA_TOC("relax");
565 }
566 }
567 }
568
569 template <class B, template <class> class C, template <class> class R>
570 friend std::ostream& operator<<(std::ostream& os, const AMG<B, C, R>& a);
571};
572
573/*---------------------------------------------------------------------------*/
574/*---------------------------------------------------------------------------*/
575
577template <class B, template <class> class C, template <class> class R>
578std::ostream& operator<<(std::ostream& os, const AMG<B, C, R>& a)
579{
580 typedef typename AMG<B, C, R>::AMGLevel level;
582
583 size_t sum_dof = 0;
584 size_t sum_nnz = 0;
585 size_t sum_mem = 0;
586
587 for (const level& lvl : a.levels) {
588 sum_dof += lvl.rows();
589 sum_nnz += lvl.nonzeros();
590 sum_mem += lvl.bytes();
591 }
592
593 os << "Number of levels: " << a.levels.size()
594 << "\nOperator complexity: " << std::fixed << std::setprecision(2)
595 << 1.0 * sum_nnz / a.levels.front().nonzeros()
596 << "\nGrid complexity: " << std::fixed << std::setprecision(2)
597 << 1.0 * sum_dof / a.levels.front().rows()
598 << "\nMemory footprint: " << human_readable_memory(sum_mem)
599 << "\n\n"
600 "level unknowns nonzeros memory\n"
601 "---------------------------------------------\n";
602
603 size_t depth = 0;
604 for (const level& lvl : a.levels) {
605 os << std::setw(5) << depth++
606 << std::setw(13) << lvl.rows()
607 << std::setw(15) << lvl.nonzeros()
608 << std::setw(12) << human_readable_memory(lvl.bytes())
609 << " (" << std::setw(5) << std::fixed << std::setprecision(2)
610 << 100.0 * lvl.nonzeros() / sum_nnz
611 << "%)" << std::endl;
612 }
613
614 return os;
615}
616
617/*---------------------------------------------------------------------------*/
618/*---------------------------------------------------------------------------*/
619
620} // namespace Arcane::Alina
621
622/*---------------------------------------------------------------------------*/
623/*---------------------------------------------------------------------------*/
624
625#endif
AMG(std::shared_ptr< build_matrix > A, const params &p=params(), const backend_params &bprm=backend_params())
Builds the AMG hierarchy for the system matrix.
Definition AMG.h:212
void cycle(const Vec1 &rhs, Vec2 &&x) const
Performs single V-cycle for the given right-hand side and solution.
Definition AMG.h:264
std::shared_ptr< matrix > system_matrix_ptr() const
Returns the system matrix from the finest level.
Definition AMG.h:291
void apply(const Vec1 &rhs, Vec2 &&x) const
Performs single V-cycle after clearing x.
Definition AMG.h:278
friend std::ostream & operator<<(std::ostream &os, const AMG< B, C, R > &a)
Sends information about the AMG hierarchy to output stream.
Definition AMG.h:578
void rebuild(const Matrix &M, const backend_params &bprm=backend_params())
Rebuild the hierarchy using the new system matrix.
Definition AMG.h:227
void rebuild(std::shared_ptr< build_matrix > A, const backend_params &bprm=backend_params())
Rebuild the hierarchy using the new system matrix.
Definition AMG.h:241
AMG(const Matrix &M, const params &p=params(), const backend_params &bprm=backend_params())
Builds the AMG hierarchy for the system matrix.
Definition AMG.h:190
Save ostream flags in constructor, restore in destructor.
Matrix class, to be used by user.
std::int32_t Int32
Type entier signé sur 32 bits.
Parameters of the method.
Definition AMG.h:99
bool allow_rebuild
Keep matrices in internal format to allow for quick rebuild of the hierarchy.
Definition AMG.h:145
Int32 npost
Number of post-relaxations.
Definition AMG.h:136
coarsening_params coarsening
Coarsening parameters.
Definition AMG.h:103
Int32 coarse_enough
Specifies when level is coarse enough to be solved directly.
Definition AMG.h:113
bool direct_coarse
Use direct solver at the coarsest level.
Definition AMG.h:121
Int32 pre_cycles
Number of cycles to make as part of preconditioning.
Definition AMG.h:142
relax_params relax
Relaxation parameters.
Definition AMG.h:104
Int32 max_levels
Maximum number of levels.
Definition AMG.h:130
Int32 ncycle
Number of cycles (1 for V-cycle, 2 for W-cycle, etc.).
Definition AMG.h:139
Int32 npre
Number of pre-relaxations.
Definition AMG.h:133