Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
DistributedBasicSolver.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/*---------------------------------------------------------------------------*/
9/*
10 * This file is based on the work on AMGCL library (version march 2026)
11 * which can be found at https://github.com/ddemidov/amgcl.
12 *
13 * Copyright (c) 2012-2022 Denis Demidov <dennis.demidov@gmail.com>
14 * SPDX-License-Identifier: MIT
15 */
16/*---------------------------------------------------------------------------*/
17/*---------------------------------------------------------------------------*/
18
19#include <iostream>
20#include <vector>
21#include <string>
22
23#include <boost/program_options.hpp>
24
25#include "arccore/alina/BuiltinBackend.h"
26#include "arccore/alina/StaticMatrix.h"
27#include "arccore/alina/Adapters.h"
28#include "arccore/alina/MessagePassingUtils.h"
29#include "arccore/alina/DistributedPreconditionedSolver.h"
30#include "arccore/alina/DistributedPreconditioner.h"
31#include "arccore/alina/DistributedSolverRuntime.h"
32
33#include "arccore/alina/IO.h"
34#include "arccore/alina/Profiler.h"
35
36#include "arcane/utils/Exception.h"
37#include "arcane/utils/PlatformUtils.h"
38#include "arcane/launcher/ArcaneLauncher.h"
39#include "arcane/core/ISubDomain.h"
40#include "arcane/utils/ITraceMng.h"
41#include "arcane/utils/IProfilingService.h"
42
43#include "AlinaSamplesCommon.h"
44
45using namespace Arcane;
46
47namespace math = Alina::math;
48
49//---------------------------------------------------------------------------
50ptrdiff_t
51assemble_poisson3d(Alina::mpi_communicator comm,
52 ptrdiff_t n, int block_size,
53 std::vector<ptrdiff_t>& ptr,
54 std::vector<ptrdiff_t>& col,
55 std::vector<double>& val,
56 std::vector<double>& rhs)
57{
58 ptrdiff_t n3 = n * n * n;
59
60 ptrdiff_t chunk = (n3 + comm.size - 1) / comm.size;
61 if (chunk % block_size != 0) {
62 chunk += block_size - chunk % block_size;
63 }
64 ptrdiff_t row_beg = std::min(n3, chunk * comm.rank);
65 ptrdiff_t row_end = std::min(n3, row_beg + chunk);
66 chunk = row_end - row_beg;
67
68 ptr.clear();
69 ptr.reserve(chunk + 1);
70 col.clear();
71 col.reserve(chunk * 7);
72 val.clear();
73 val.reserve(chunk * 7);
74
75 rhs.resize(chunk);
76 std::fill(rhs.begin(), rhs.end(), 1.0);
77
78 const double h2i = (n - 1) * (n - 1);
79 ptr.push_back(0);
80
81 for (ptrdiff_t idx = row_beg; idx < row_end; ++idx) {
82 ptrdiff_t k = idx / (n * n);
83 ptrdiff_t j = (idx / n) % n;
84 ptrdiff_t i = idx % n;
85
86 if (k > 0) {
87 col.push_back(idx - n * n);
88 val.push_back(-h2i);
89 }
90
91 if (j > 0) {
92 col.push_back(idx - n);
93 val.push_back(-h2i);
94 }
95
96 if (i > 0) {
97 col.push_back(idx - 1);
98 val.push_back(-h2i);
99 }
100
101 col.push_back(idx);
102 val.push_back(6 * h2i);
103
104 if (i + 1 < n) {
105 col.push_back(idx + 1);
106 val.push_back(-h2i);
107 }
108
109 if (j + 1 < n) {
110 col.push_back(idx + n);
111 val.push_back(-h2i);
112 }
113
114 if (k + 1 < n) {
115 col.push_back(idx + n * n);
116 val.push_back(-h2i);
117 }
118
119 ptr.push_back(col.size());
120 }
121
122 return chunk;
123}
124
125//---------------------------------------------------------------------------
126template <class Backend, class Matrix>
127std::shared_ptr<Alina::DistributedMatrix<Backend>>
128partition(Alina::mpi_communicator comm, const Matrix& Astrip,
129 typename Backend::vector& rhs, const typename Backend::params& bprm,
130 Alina::eMatrixPartitionerType ptype, int block_size = 1)
131{
132 auto& prof = Alina::Profiler::globalProfiler();
133 typedef typename Backend::value_type val_type;
134 typedef typename Alina::math::rhs_of<val_type>::type rhs_type;
135 typedef Alina::DistributedMatrix<Backend> DMatrix;
136
137 auto A = std::make_shared<DMatrix>(comm, Astrip);
138
139 if (comm.size == 1 || ptype == Alina::eMatrixPartitionerType::merge)
140 return A;
141
142 prof.tic("partition");
144 prm.put("type", ptype);
146
147 auto I = part(*A, block_size);
148 auto J = transpose(*I);
149 A = product(*J, *product(*A, *I));
150
151 Alina::numa_vector<rhs_type> new_rhs(J->loc_rows());
152
153 J->move_to_backend(bprm);
154
155 Alina::backend::spmv(1, *J, rhs, 0, new_rhs);
156 rhs.swap(new_rhs);
157 prof.toc("partition");
158
159 return A;
160}
161
162//---------------------------------------------------------------------------
163void solve_scalar(Alina::mpi_communicator comm,
164 ptrdiff_t chunk,
165 const std::vector<ptrdiff_t>& ptr,
166 const std::vector<ptrdiff_t>& col,
167 const std::vector<double>& val,
168 const Alina::PropertyTree& prm,
169 const std::vector<double>& f,
170 Alina::eMatrixPartitionerType ptype)
171{
172 auto& prof = Alina::Profiler::globalProfiler();
173 //using Backend = Alina::BuiltinBackend<double>;
174
175 using BackendValueType = double;
177
178 std::cout << "Using scalar solve ptr_size=" << sizeof(ptrdiff_t)
179 << " ptr_type_size=" << sizeof(Backend::ptr_type)
180 << " col_type_size=" << sizeof(Backend::col_type)
181 << " value_type_size=" << sizeof(Backend::value_type)
182 << "\n";
183
184 typedef Alina::DistributedMatrix<Backend> DMatrix;
185
187 using RelaxationType = Alina::DistributedSPAI0Relaxation<Backend>;
188 // Si on veut tester les backend dynamiques:
189 //using CoarseningType = Alina::DistributedCoarseningRuntime<Backend>;
190 //using RelaxationType = Alina::DistributedRelaxationRuntime<Backend>,
191
192 using AMGPrecondType = Alina::DistributedAMG<Backend, CoarseningType, RelaxationType,
195
197
198 typename Backend::params bprm;
199
201
202 auto get_distributed_matrix = [&]() {
203 auto t = prof.scoped_tic("distributed matrix");
204 std::shared_ptr<DMatrix> A;
205
206 if (ptype != Alina::eMatrixPartitionerType::merge) {
207 A = partition<Backend>(comm,
208 std::tie(chunk, ptr, col, val), rhs, bprm, ptype,
209 prm.get("precond.coarsening.aggr.block_size", 1));
210 chunk = A->loc_rows();
211 }
212 else {
213 A = std::make_shared<DMatrix>(comm, std::tie(chunk, ptr, col, val));
214 }
215
216 return A;
217 };
218
219 std::shared_ptr<DMatrix> A;
220 std::shared_ptr<Solver> solve;
221
222 {
223 auto t = prof.scoped_tic("setup");
224 A = get_distributed_matrix();
225 solve = std::make_shared<Solver>(comm, A, prm, bprm);
227 solve->prm.get(prm2);
228 std::cout << "SOLVER parameters=" << prm2 << "\n";
229 }
230
231 if (comm.rank == 0) {
232 std::cout << "SolverInfo:\n";
233 std::cout << *solve << std::endl;
234 }
235
236 if (prm.get("precond.allow_rebuild", false)) {
237 if (comm.rank == 0) {
238 std::cout << "Rebuilding the preconditioner..." << std::endl;
239 }
240
241 {
242 auto t = prof.scoped_tic("rebuild");
243 A = get_distributed_matrix();
244 solve->precond().rebuild(A);
245 }
246
247 if (comm.rank == 0) {
248 std::cout << *solve << std::endl;
249 }
250 }
251
253
254 prof.tic("solve");
255 Alina::SolverResult r = (*solve)(rhs, x);
256 prof.toc("solve");
257
258 if (comm.rank == 0) {
259 std::cout << "Iterations: " << r.nbIteration() << std::endl
260 << "Error: " << r.residual() << std::endl
261 << prof << std::endl;
262 }
263}
264
265//---------------------------------------------------------------------------
266int main2(const Alina::SampleMainContext& ctx, int argc, char* argv[])
267{
268 ITraceMng* tm = ctx.traceMng();
269 auto& prof = Alina::Profiler::globalProfiler();
270
271 //Alina::mpi_init_thread mpi(&argc, &argv);
272 Alina::mpi_communicator comm(MPI_COMM_WORLD);
273
274 tm->info() << "World size: " << comm.size;
275
276 // Read configuration from command line
277 namespace po = boost::program_options;
278 po::options_description desc("Options");
279
280 auto default_partitioner_type = Alina::eMatrixPartitionerType::merge;
281#if defined(ARCCORE_ALINA_HAVE_PARMETIS)
282 default_partitioner_type = Alina::eMatrixPartitionerType::parmetis;
283#endif
284
285 desc.add_options()("help,h", "show help")("matrix,A",
286 po::value<std::string>(),
287 "System matrix in the MatrixMarket format. "
288 "When not specified, a Poisson problem in 3D unit cube is assembled. ");
289 desc.add_options()("partitioner,r",
290 po::value<Alina::eMatrixPartitionerType>()->default_value(
291 default_partitioner_type),
292 "Repartition the system matrix");
293 desc.add_options()("size,n",
294 po::value<ptrdiff_t>()->default_value(32),
295 "domain size");
296 desc.add_options()("prm-file,P", po::value<std::string>(),
297 "Parameter file in json format. ");
298 desc.add_options()("prm,p",
299 po::value<std::vector<std::string>>()->multitoken(),
300 "Parameters specified as name=value pairs. "
301 "May be provided multiple times. Examples:\n"
302 " -p solver.tol=1e-3\n"
303 " -p precond.coarse_enough=300");
304 desc.add_options()("test-rebuild",
305 po::bool_switch()->default_value(false),
306 "When specified, try to rebuild the solver before solving. ");
307
308 po::positional_options_description p;
309 p.add("prm", -1);
310
311 po::variables_map vm;
312 po::store(po::command_line_parser(argc, argv).options(desc).positional(p).run(), vm);
313 po::notify(vm);
314
315 if (vm.count("help")) {
316 if (comm.rank == 0)
317 std::cout << desc << std::endl;
318 return 0;
319 }
320
322 if (vm.count("prm-file")) {
323 prm.read_json(vm["prm-file"].as<std::string>());
324 }
325
326 if (vm.count("prm")) {
327 for (const std::string& v : vm["prm"].as<std::vector<std::string>>()) {
328 tm->info() << "PUT_KEY_VALUE v=" << v;
329 prm.putKeyValue(v);
330 }
331 }
332
333 ptrdiff_t n;
334 std::vector<ptrdiff_t> ptr;
335 std::vector<ptrdiff_t> col;
336 std::vector<double> val;
337 std::vector<double> rhs;
338
339 Alina::eMatrixPartitionerType ptype = vm["partitioner"].as<Alina::eMatrixPartitionerType>();
340
341 prof.tic("assemble");
342 Int64 matrix_size = vm["size"].as<ptrdiff_t>();
343 tm->info() << "Matrix size=" << matrix_size;
344 n = assemble_poisson3d(comm, matrix_size, 1, ptr, col, val, rhs);
345 prof.toc("assemble");
346
347 if (vm["test-rebuild"].as<bool>()) {
348 prm.put("precond.allow_rebuild", true);
349 }
350
351 solve_scalar(comm, n, ptr, col, val, prm, rhs, ptype);
352 return 0;
353}
354
355/*---------------------------------------------------------------------------*/
356/*---------------------------------------------------------------------------*/
357
358int main(int argc, char* argv[])
359{
360 return Arcane::Alina::SampleMainContext::execMain(main2, argc, argv);
361}
362
363/*---------------------------------------------------------------------------*/
364/*---------------------------------------------------------------------------*/
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
Interface du gestionnaire de traces.
virtual TraceMessage info()=0
Flot pour un message d'information.
Matrix class, to be used by user.
Espace de nom pour les fonctions mathématiques.
Definition MathUtils.h:41
__host__ __device__ Real3x3 transpose(const Real3x3 &t)
Transpose la matrice.
Definition MathUtils.h:258
-*- 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.