Arcane  4.1.11.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
TestDistributedSolver.cc
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#include <gtest/gtest.h>
9
10#include "arccore/alina/BuiltinBackend.h"
11#include "arccore/alina/StaticMatrix.h"
12#include "arccore/alina/Adapters.h"
13#include "arccore/alina/MessagePassingUtils.h"
14#include "arccore/alina/DistributedPreconditionedSolver.h"
15#include "arccore/alina/DistributedPreconditioner.h"
16#include "arccore/alina/DistributedSolverRuntime.h"
17
18#include "TestMainMpi.h"
19
20using namespace Arcane;
21
22namespace math = Alina::math;
23
24//---------------------------------------------------------------------------
25ptrdiff_t
26assemble_poisson3d(Alina::mpi_communicator comm,
27 ptrdiff_t n, int block_size,
28 std::vector<ptrdiff_t>& ptr,
29 std::vector<ptrdiff_t>& col,
30 std::vector<double>& val,
31 std::vector<double>& rhs)
32{
33 ptrdiff_t n3 = n * n * n;
34
35 ptrdiff_t chunk = (n3 + comm.size - 1) / comm.size;
36 if (chunk % block_size != 0) {
37 chunk += block_size - chunk % block_size;
38 }
39 ptrdiff_t row_beg = std::min(n3, chunk * comm.rank);
40 ptrdiff_t row_end = std::min(n3, row_beg + chunk);
41 chunk = row_end - row_beg;
42
43 ptr.clear();
44 ptr.reserve(chunk + 1);
45 col.clear();
46 col.reserve(chunk * 7);
47 val.clear();
48 val.reserve(chunk * 7);
49
50 rhs.resize(chunk);
51 std::fill(rhs.begin(), rhs.end(), 1.0);
52
53 const double h2i = (n - 1) * (n - 1);
54 ptr.push_back(0);
55
56 for (ptrdiff_t idx = row_beg; idx < row_end; ++idx) {
57 ptrdiff_t k = idx / (n * n);
58 ptrdiff_t j = (idx / n) % n;
59 ptrdiff_t i = idx % n;
60
61 if (k > 0) {
62 col.push_back(idx - n * n);
63 val.push_back(-h2i);
64 }
65
66 if (j > 0) {
67 col.push_back(idx - n);
68 val.push_back(-h2i);
69 }
70
71 if (i > 0) {
72 col.push_back(idx - 1);
73 val.push_back(-h2i);
74 }
75
76 col.push_back(idx);
77 val.push_back(6 * h2i);
78
79 if (i + 1 < n) {
80 col.push_back(idx + 1);
81 val.push_back(-h2i);
82 }
83
84 if (j + 1 < n) {
85 col.push_back(idx + n);
86 val.push_back(-h2i);
87 }
88
89 if (k + 1 < n) {
90 col.push_back(idx + n * n);
91 val.push_back(-h2i);
92 }
93
94 ptr.push_back(col.size());
95 }
96
97 return chunk;
98}
99
100//---------------------------------------------------------------------------
101
102void solve_scalar(Alina::mpi_communicator comm,
103 ptrdiff_t chunk,
104 const std::vector<ptrdiff_t>& ptr,
105 const std::vector<ptrdiff_t>& col,
106 const std::vector<double>& val,
107 const Alina::PropertyTree& prm,
108 const std::vector<double>& f)
109{
110 auto& prof = Alina::Profiler::globalProfiler();
111 //using Backend = Alina::BuiltinBackend<double>;
112
113 using BackendValueType = double;
115
116 std::cout << "Using scalar solve ptr_size=" << sizeof(ptrdiff_t)
117 << " ptr_type_size=" << sizeof(Backend::ptr_type)
118 << " col_type_size=" << sizeof(Backend::col_type)
119 << " value_type_size=" << sizeof(Backend::value_type)
120 << "\n";
121
122 typedef Alina::DistributedMatrix<Backend> DMatrix;
123
125 using RelaxationType = Alina::DistributedSPAI0Relaxation<Backend>;
126 // Si on veut tester les backend dynamiques:
127 //using CoarseningType = Alina::DistributedCoarseningRuntime<Backend>;
128 //using RelaxationType = Alina::DistributedRelaxationRuntime<Backend>,
129
130 using AMGPrecondType = Alina::DistributedAMG<Backend, CoarseningType, RelaxationType,
133
135
136 typename Backend::params bprm;
137
139
140 std::shared_ptr<DMatrix> A;
141 std::shared_ptr<Solver> solve;
142
143 {
144 auto t = prof.scoped_tic("setup");
145 A = std::make_shared<DMatrix>(comm, std::tie(chunk, ptr, col, val));
146 solve = std::make_shared<Solver>(comm, A, prm, bprm);
148 solve->prm.get(prm2);
149 std::cout << "SOLVER parameters=" << prm2 << "\n";
150 }
151
152 if (comm.rank == 0) {
153 std::cout << "SolverInfo:\n";
154 std::cout << *solve << std::endl;
155 }
156
158
159 prof.tic("solve");
160 Alina::SolverResult r = (*solve)(rhs, x);
161 prof.toc("solve");
162
163 if (comm.rank == 0) {
164 std::cout << "Iterations: " << r.nbIteration() << std::endl
165 << "Error: " << r.residual() << std::endl
166 << prof << std::endl;
167 }
168}
169
170/*---------------------------------------------------------------------------*/
171/*---------------------------------------------------------------------------*/
172
173TEST(alina_test_mpi, BasicSolver)
174{
175 Alina::mpi_communicator comm(AlinaTest::global_mpi_comm_world);
176
177 std::cout << "World size: " << comm.size << "\n";
178
180
181 ptrdiff_t n;
182 std::vector<ptrdiff_t> ptr;
183 std::vector<ptrdiff_t> col;
184 std::vector<double> val;
185 std::vector<double> rhs;
186
187 Int64 matrix_size = 32;
188 std::cout << "Matrix size=" << matrix_size << "\n";
189 n = assemble_poisson3d(comm, matrix_size, 1, ptr, col, val, rhs);
190
191 solve_scalar(comm, n, ptr, col, val, prm, rhs);
192}
193
194/*---------------------------------------------------------------------------*/
195/*---------------------------------------------------------------------------*/
Runtime wrapper for distributed direct solvers.
Distributed Matrix using message passing.
Iterative solver wrapper for distributed linear systems.
Result of a solving.
Definition AlinaUtils.h:52
NUMA-aware vector container.
Definition NumaVector.h:42
Espace de nom pour les fonctions mathématiques.
Definition MathUtils.h:41
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Alina::detail::empty_params params
Distributed smoothed aggregation coarsening scheme.
Runtime-configurable wrapper around matrix partitioner.
Convenience wrapper around MPI_Comm.