Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
DistributedAMG.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/* DistributedAMG.h (C) 2000-2026 */
9/* */
10/* Distributed memory AMG preconditioner. */
11/*---------------------------------------------------------------------------*/
12#ifndef ARCCORE_ALINA_DISTRIBUTEDAMG_H
13#define ARCCORE_ALINA_DISTRIBUTEDAMG_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 <iostream>
27#include <iomanip>
28#include <list>
29#include <memory>
30
31#include "arccore/alina/BackendInterface.h"
32#include "arccore/alina/MessagePassingUtils.h"
33#include "arccore/alina/DistributedMatrix.h"
34#include "arccore/alina/DistributedSkylineLUDirectSolver.h"
35#include "arccore/alina/SimpleMatrixPartitioner.h"
36
37/*---------------------------------------------------------------------------*/
38/*---------------------------------------------------------------------------*/
39
40namespace Arcane::Alina
41{
42
43/*---------------------------------------------------------------------------*/
44/*---------------------------------------------------------------------------*/
45
46template <class Backend,
47 class Coarsening,
48 class Relaxation,
50 class Repartition = SimpleMatrixPartitioner<Backend>>
51class DistributedAMG
52{
53 public:
54
55 using backend_type = Backend;
56 using BackendType = Backend;
57
58 typedef typename Backend::params backend_params;
59 typedef typename Backend::value_type value_type;
60 typedef typename math::scalar_of<value_type>::type scalar_type;
61 typedef DistributedMatrix<Backend> matrix;
62 typedef typename Backend::vector vector;
63
64 struct params
65 {
66 typedef typename Coarsening::params coarsening_params;
67 typedef typename Relaxation::params relax_params;
68 typedef typename DirectSolver::params direct_params;
69 typedef typename Repartition::params repart_params;
70
71 coarsening_params coarsening;
72 relax_params relax;
73 direct_params direct;
74 repart_params repart;
75
83 Int32 coarse_enough = DirectSolver::coarse_enough();
84
91 bool direct_coarse = true;
92
100 Int32 max_levels = std::numeric_limits<Int32>::max();
101
104
107
110
113
115 bool allow_rebuild = false;
116
117 params() = default;
118
119 params(const PropertyTree& p)
120 : ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p, coarsening)
121 , ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p, relax)
122 , ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p, direct)
123 , ARCCORE_ALINA_PARAMS_IMPORT_CHILD(p, repart)
124 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, coarse_enough)
125 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, direct_coarse)
126 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, max_levels)
127 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, npre)
128 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, npost)
129 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, ncycle)
130 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, pre_cycles)
131 , ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, allow_rebuild)
132 {
133 p.check_params({ "coarsening", "relax", "direct", "repart", "coarse_enough", "direct_coarse", "max_levels", "npre", "npost", "ncycle", "pre_cycles", "allow_rebuild" });
134
135 Alina::precondition(max_levels > 0, "max_levels should be positive");
136 }
137
138 void get(Alina::PropertyTree& p, const std::string& path = "") const
139 {
140 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, coarsening);
141 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, relax);
142 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, direct);
143 ARCCORE_ALINA_PARAMS_EXPORT_CHILD(p, path, repart);
144 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, coarse_enough);
145 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, direct_coarse);
146 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, max_levels);
147 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, npre);
148 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, npost);
149 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, ncycle);
150 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, pre_cycles);
151 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, allow_rebuild);
152 }
153 } prm;
154
155 template <class Matrix>
156 DistributedAMG(mpi_communicator comm,
157 const Matrix& A,
158 const params& prm = params(),
159 const backend_params& bprm = backend_params())
160 : prm(prm)
161 , comm(comm)
162 , repart(prm.repart)
163 {
164 init(std::make_shared<matrix>(comm, A, backend::nbRow(A)), bprm);
165 }
166
167 DistributedAMG(mpi_communicator comm,
168 std::shared_ptr<matrix> A,
169 const params& prm = params(),
170 const backend_params& bprm = backend_params())
171 : prm(prm)
172 , comm(comm)
173 , repart(prm.repart)
174 {
175 init(A, bprm);
176 }
177
184 template <class Matrix>
185 void rebuild(const Matrix& M,
186 const backend_params& bprm = backend_params())
187 {
188 rebuild(std::make_shared<matrix>(comm, M, backend::nbRow(M)), bprm);
189 }
190
191 template <class OtherBackend>
192 typename std::enable_if<!std::is_same<Backend, OtherBackend>::value, void>::type
193 rebuild(std::shared_ptr<DistributedMatrix<OtherBackend>> A,
194 const backend_params& bprm = backend_params())
195 {
196 return rebuild(std::make_shared<matrix>(*A), bprm);
197 }
198
199 void rebuild(std::shared_ptr<matrix> A,
200 const backend_params& bprm = backend_params())
201 {
202 precondition(prm.allow_rebuild, "allow_rebuild is not set!");
203 precondition(A->glob_rows() == system_matrix().glob_rows() &&
204 A->glob_cols() == system_matrix().glob_cols(),
205 "Matrix dimensions differ from the original ones!");
206
207 ARCCORE_ALINA_TIC("rebuild");
208 this->A = A;
209 Coarsening C(prm.coarsening);
210 for (auto& level : levels) {
211 A = level.rebuild(A, C, prm, bprm);
212 }
213 this->A->move_to_backend(bprm);
214 ARCCORE_ALINA_TOC("rebuild");
215 }
216
217 template <class Vec1, class Vec2>
218 void cycle(const Vec1& rhs, Vec2&& x) const
219 {
220 cycle(levels.begin(), rhs, x);
221 }
222
223 template <class Vec1, class Vec2>
224 void apply(const Vec1& rhs, Vec2&& x) const
225 {
226 if (prm.pre_cycles) {
227 backend::clear(x);
228 for (unsigned i = 0; i < prm.pre_cycles; ++i)
229 cycle(levels.begin(), rhs, x);
230 }
231 else {
232 backend::copy(rhs, x);
233 }
234 }
235
237 std::shared_ptr<matrix> system_matrix_ptr() const
238 {
239 return A;
240 }
241
242 const matrix& system_matrix() const
243 {
244 return *system_matrix_ptr();
245 }
246
247 private:
248
249 struct DistributedAMGLevel
250 {
251 ptrdiff_t nrows, nnz;
252 int active_procs;
253
254 std::shared_ptr<matrix> A, P, R;
255 std::shared_ptr<vector> f, u, t;
256 std::shared_ptr<Relaxation> relax;
257 std::shared_ptr<DirectSolver> solve;
258
259 DistributedAMGLevel() = default;
260
261 DistributedAMGLevel(std::shared_ptr<matrix> a,
262 params& prm,
263 const backend_params& bprm,
264 bool direct = false)
265 : nrows(a->glob_rows())
266 , nnz(a->glob_nonzeros())
267 , f(Backend::create_vector(a->loc_rows(), bprm))
268 , u(Backend::create_vector(a->loc_rows(), bprm))
269 {
270 int active = (a->loc_rows() > 0);
271 active_procs = a->comm().reduceSum(active);
272
273 sort_rows(*a);
274
275 if (direct) {
276 ARCCORE_ALINA_TIC("direct solver");
277 solve = std::make_shared<DirectSolver>(a->comm(), *a, prm.direct);
278 ARCCORE_ALINA_TOC("direct solver");
279 }
280 else {
281 A = a;
282 t = Backend::create_vector(a->loc_rows(), bprm);
283
284 ARCCORE_ALINA_TIC("relaxation");
285 relax = std::make_shared<Relaxation>(*a, prm.relax, bprm);
286 ARCCORE_ALINA_TOC("relaxation");
287 }
288 }
289
290 std::shared_ptr<matrix> step_down(Coarsening& C, const Repartition& repart)
291 {
292 ARCCORE_ALINA_TIC("transfer operators");
293 std::tie(P, R) = C.transfer_operators(*A);
294
295 ARCCORE_ALINA_TIC("sort");
296 sort_rows(*P);
297 sort_rows(*R);
298 ARCCORE_ALINA_TOC("sort");
299
300 ARCCORE_ALINA_TOC("transfer operators");
301
302 if (P->glob_cols() == 0) {
303 // Zero-sized coarse level in AMG (diagonal matrix?)
304 return std::shared_ptr<matrix>();
305 }
306
307 ARCCORE_ALINA_TIC("coarse operator");
308 auto Ac = C.coarse_operator(*A, *P, *R);
309 ARCCORE_ALINA_TOC("coarse operator");
310
311 if (repart.is_needed(*Ac)) {
312 ARCCORE_ALINA_TIC("partition");
313 auto I = repart(*Ac, block_size(C));
314 auto J = transpose(*I);
315
316 P = product(*P, *I);
317 R = product(*J, *R);
318 Ac = product(*J, *product(*Ac, *I));
319 ARCCORE_ALINA_TOC("partition");
320 }
321
322 return Ac;
323 }
324
325 std::shared_ptr<matrix> rebuild(std::shared_ptr<matrix> A,
326 const Coarsening& C,
327 const params& prm,
328 const backend_params& bprm)
329 {
330 if (relax) {
331 relax = std::make_shared<Relaxation>(*A, prm.relax, bprm);
332 std::cout << "DistributedAMG: relaxation=" << relax << "\n";
333 }
334
335 if (solve) {
336 solve = std::make_shared<DirectSolver>(A->comm(), *A, prm.direct);
337 }
338
339 if (this->A) {
340 this->A = A;
341 }
342
343 if (P && R) {
344 A = C.coarse_operator(*A, *P, *R);
345 }
346
347 if (this->A) {
348 this->A->move_to_backend(bprm);
349 }
350
351 return A;
352 }
353
354 void move_to_backend(const backend_params& bprm, bool keep_src = false)
355 {
356 ARCCORE_ALINA_TIC("move to backend");
357 if (A)
358 A->move_to_backend(bprm);
359 if (P)
360 P->move_to_backend(bprm, keep_src);
361 if (R)
362 R->move_to_backend(bprm, keep_src);
363 ARCCORE_ALINA_TOC("move to backend");
364 }
365
366 ptrdiff_t rows() const
367 {
368 return nrows;
369 }
370
371 ptrdiff_t nonzeros() const
372 {
373 return nnz;
374 }
375 };
376
377 typedef typename std::list<DistributedAMGLevel>::const_iterator level_iterator;
378
379 mpi_communicator comm;
380 std::shared_ptr<matrix> A;
381 Repartition repart;
382 std::list<DistributedAMGLevel> levels;
383
384 void init(std::shared_ptr<matrix> A, const backend_params& bprm)
385 {
386 A->comm().check(A->glob_rows() == A->glob_cols(), "Matrix should be square!");
387
388 this->A = A;
389 Coarsening C(prm.coarsening);
390 //std::cout << "DistributedAMGInit: Coarsening=" << C << "\n";
391 bool need_coarse = true;
392
393 while (A->glob_rows() > prm.coarse_enough) {
394 levels.push_back(DistributedAMGLevel(A, prm, bprm));
395
396 if (levels.size() >= prm.max_levels) {
397 levels.back().move_to_backend(bprm, prm.allow_rebuild);
398 break;
399 }
400
401 A = levels.back().step_down(C, repart);
402 levels.back().move_to_backend(bprm, prm.allow_rebuild);
403
404 if (!A) {
405 // Zero-sized coarse level. Probably the system matrix on
406 // this level is diagonal, should be easily solvable with a
407 // couple of smoother iterations.
408 need_coarse = false;
409 break;
410 }
411 }
412
413 if (!A || A->glob_rows() > prm.coarse_enough) {
414 // The coarse matrix is still too big to be solved directly.
415 need_coarse = false;
416 }
417
418 if (A && need_coarse) {
419 levels.push_back(DistributedAMGLevel(A, prm, bprm, prm.direct_coarse));
420 levels.back().move_to_backend(bprm, prm.allow_rebuild);
421 }
422
423 ARCCORE_ALINA_TIC("move to backend");
424 this->A->move_to_backend(bprm, prm.allow_rebuild);
425 ARCCORE_ALINA_TOC("move to backend");
426 }
427
428 template <class Vec1, class Vec2>
429 void cycle(level_iterator lvl, const Vec1& rhs, Vec2& x) const
430 {
431 level_iterator nxt = lvl, end = levels.end();
432 ++nxt;
433
434 if (nxt == end) {
435 if (lvl->solve) {
436 ARCCORE_ALINA_TIC("direct solver");
437 (*lvl->solve)(rhs, x);
438 ARCCORE_ALINA_TOC("direct solver");
439 }
440 else {
441 ARCCORE_ALINA_TIC("relax");
442 for (size_t i = 0; i < prm.npre; ++i)
443 lvl->relax->apply_pre(*lvl->A, rhs, x, *lvl->t);
444 for (size_t i = 0; i < prm.npost; ++i)
445 lvl->relax->apply_post(*lvl->A, rhs, x, *lvl->t);
446 ARCCORE_ALINA_TOC("relax");
447 }
448 }
449 else {
450 for (size_t j = 0; j < prm.ncycle; ++j) {
451 ARCCORE_ALINA_TIC("relax");
452 for (size_t i = 0; i < prm.npre; ++i)
453 lvl->relax->apply_pre(*lvl->A, rhs, x, *lvl->t);
454 ARCCORE_ALINA_TOC("relax");
455
456 backend::residual(rhs, *lvl->A, x, *lvl->t);
457
458 backend::spmv(math::identity<scalar_type>(), *lvl->R, *lvl->t, math::zero<scalar_type>(), *nxt->f);
459
460 backend::clear(*nxt->u);
461 cycle(nxt, *nxt->f, *nxt->u);
462
463 backend::spmv(math::identity<scalar_type>(), *lvl->P, *nxt->u, math::identity<scalar_type>(), x);
464
465 ARCCORE_ALINA_TIC("relax");
466 for (size_t i = 0; i < prm.npost; ++i)
467 lvl->relax->apply_post(*lvl->A, rhs, x, *lvl->t);
468 ARCCORE_ALINA_TOC("relax");
469 }
470 }
471 }
472
473 template <class B, class C, class R, class D, class I>
474 friend std::ostream& operator<<(std::ostream& os, const DistributedAMG<B, C, R, D, I>& a);
475};
476
477template <class B, class C, class R, class D, class I>
478std::ostream& operator<<(std::ostream& os, const DistributedAMG<B, C, R, D, I>& a)
479{
482
483 size_t sum_dof = 0;
484 size_t sum_nnz = 0;
485
486 for (const level& lvl : a.levels) {
487 sum_dof += lvl.rows();
488 sum_nnz += lvl.nonzeros();
489 }
490
491 os << "Preconditioner: DistributedAMG\n";
492 //os << "Coarsening: " << a.prm.coarsening.type << "\n";
493 //os << "Relaxation: " << a.prm.relax.type << "\n";
494
495 os << "Number of levels: " << a.levels.size()
496 << "\nOperator complexity: " << std::fixed << std::setprecision(2)
497 << 1.0 * sum_nnz / a.levels.front().nonzeros()
498 << "\nGrid complexity: " << std::fixed << std::setprecision(2)
499 << 1.0 * sum_dof / a.levels.front().rows()
500 << "\n\nlevel unknowns nonzeros\n"
501 << "---------------------------------\n";
502
503 size_t depth = 0;
504 for (const level& lvl : a.levels) {
505 os << std::setw(5) << depth++
506 << std::setw(13) << lvl.rows()
507 << std::setw(15) << lvl.nonzeros() << " ("
508 << std::setw(5) << std::fixed << std::setprecision(2)
509 << 100.0 * lvl.nonzeros() / sum_nnz
510 << "%) [" << lvl.active_procs << "]" << std::endl;
511 }
512
513 return os;
514}
515
516/*---------------------------------------------------------------------------*/
517/*---------------------------------------------------------------------------*/
518
519} // namespace Arcane::Alina
520
521/*---------------------------------------------------------------------------*/
522/*---------------------------------------------------------------------------*/
523
524#endif
void rebuild(const Matrix &M, const backend_params &bprm=backend_params())
Rebuild the hierarchy using the new system matrix.
std::shared_ptr< matrix > system_matrix_ptr() const
Returns the system matrix from the finest level.
Distributed Matrix using message passing.
Provides distributed direct solver interface for Skyline LU solver.
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.
repart_params repart
Repartition parameters.
Int32 coarse_enough
Specifies when level is coarse enough to be solved directly.
relax_params relax
Relaxation parameters.
bool allow_rebuild
Keep matrices in internal format to allow for quick rebuild of the hierarchy.
Int32 npre
Number of pre-relaxations.
Int32 npost
Number of post-relaxations.
Int32 ncycle
Number of cycles (1 for V-cycle, 2 for W-cycle, etc.).
Int32 max_levels
Maximum number of levels.
coarsening_params coarsening
Coarsening parameters.
direct_params direct
Direct solver parameters.
Int32 pre_cycles
Number of cycles to make as part of preconditioning.
bool direct_coarse
Use direct solver at the coarsest level.
Simple matrix partitioner merging consecutive domains together.
Convenience wrapper around MPI_Comm.