Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
Solver.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 <string>
21#include <random>
22
23// To remove warnings about deprecated Eigen usage.
24#pragma GCC diagnostic ignored "-Wdeprecated-copy"
25#pragma GCC diagnostic ignored "-Wint-in-bool-context"
26
27#include <boost/program_options.hpp>
28#include <boost/range/iterator_range.hpp>
29#include <boost/preprocessor/seq/for_each.hpp>
30
31#if defined(SOLVER_BACKEND_CUDA)
32#include "arccore/alina/CudaBackend.h"
33#include "arccore/alina/relaxation_cusparse_ilu0.h"
35#elif defined(SOLVER_BACKEND_EIGEN)
36#include "arccore/alina/EigenBackend.h"
38#else
39#ifndef SOLVER_BACKEND_BUILTIN
40#define SOLVER_BACKEND_BUILTIN
41#endif
42#include "arccore/alina/BuiltinBackend.h"
43#include "arccore/alina/StaticMatrix.h"
44#include "arccore/alina/Adapters.h"
45// Use 32 bit indexing for backend.
47//using Backend = Arcane::Alina::BuiltinBackend<double>;
48#endif
49
50#include <arcane/utils/PlatformUtils.h>
51#include <arcane/utils/String.h>
52#include <arcane/utils/Convert.h>
53
54#include "arccore/alina/RelaxationRuntime.h"
55#include "arccore/alina/CoarseningRuntime.h"
56#include "arccore/alina/SolverRuntime.h"
57#include "arccore/alina/PreconditionerRuntime.h"
58#include "arccore/alina/PreconditionedSolver.h"
59#include "arccore/alina/AMG.h"
60#include "arccore/alina/Adapters.h"
61#include "arccore/alina/IO.h"
62
63#include "arccore/alina/Profiler.h"
64
65#include "SampleProblemCommon.h"
66
67#ifndef ARCCORE_ALINA_BLOCK_SIZES
68#define ARCCORE_ALINA_BLOCK_SIZES (3)(4)
69#endif
70
71using namespace Arcane;
72
73using Alina::precondition;
74
75#ifdef SOLVER_BACKEND_BUILTIN
76extern "C++" void
77_doHypreSolver(int nb_row,
78 std::vector<ptrdiff_t> const& ptr,
79 std::vector<ptrdiff_t> const& col,
80 std::vector<double> const& val,
81 std::vector<double> const& rhs,
82 std::vector<double>& x,
83 int argc, char* argv[]);
84#endif
85
86#ifdef SOLVER_BACKEND_BUILTIN
87//---------------------------------------------------------------------------
88template <int B> Alina::SolverResult
89block_solve(const Alina::PropertyTree& prm,
90 size_t rows,
91 std::vector<ptrdiff_t> const& ptr,
92 std::vector<ptrdiff_t> const& col,
93 std::vector<double> const& val,
94 std::vector<double> const& rhs,
95 std::vector<double>& x,
96 bool reorder)
97{
98 auto& prof = Alina::Profiler::globalProfiler();
99
100 typedef Alina::StaticMatrix<double, B, B> value_type;
101 typedef Alina::StaticMatrix<double, B, 1> rhs_type;
102 typedef Alina::BuiltinBackend<value_type> BBackend;
103
105
106 auto As = std::tie(rows, ptr, col, val);
107 auto Ab = Alina::adapter::block_matrix<value_type>(As);
108
109 std::tuple<size_t, double> info;
110
111 if (reorder) {
112 prof.tic("reorder");
114 prof.toc("reorder");
115
116 prof.tic("setup");
117 Solver solve(perm(Ab), prm);
118 prof.toc("setup");
119
120 std::cout << solve << std::endl;
121
122 rhs_type const* fptr = reinterpret_cast<rhs_type const*>(&rhs[0]);
123 rhs_type* xptr = reinterpret_cast<rhs_type*>(&x[0]);
124
126 Alina::numa_vector<rhs_type> X(perm(SmallSpan<rhs_type>(xptr, rows / B)));
127
128 prof.tic("solve");
129 info = solve(F, X);
130 prof.toc("solve");
131
132 perm.inverse(X, xptr);
133 }
134 else {
135 prof.tic("setup");
136 Solver solve(Ab, prm);
137 prof.toc("setup");
138
139 std::cout << solve << std::endl;
140
141 rhs_type const* fptr = reinterpret_cast<rhs_type const*>(&rhs[0]);
142 rhs_type* xptr = reinterpret_cast<rhs_type*>(&x[0]);
143
144 Alina::numa_vector<rhs_type> F(fptr, fptr + rows / B);
145 Alina::numa_vector<rhs_type> X(xptr, xptr + rows / B);
146
147 prof.tic("solve");
148 info = solve(F, X);
149 prof.toc("solve");
150
151 std::copy(X.data(), X.data() + X.size(), xptr);
152 }
153
154 return info;
155}
156#endif
157
158//---------------------------------------------------------------------------
160scalar_solve(const Alina::PropertyTree& prm,
161 size_t rows,
162 std::vector<ptrdiff_t> const& ptr,
163 std::vector<ptrdiff_t> const& col,
164 std::vector<double> const& val,
165 std::vector<double> const& rhs,
166 std::vector<double>& x,
167 bool reorder)
168{
169 std::cout << "Using scalar solve ptr_size=" << sizeof(ptrdiff_t)
170 << " ptr_type_size=" << sizeof(Backend::ptr_type)
171 << " col_type_size=" << sizeof(Backend::col_type)
172 << " value_type_size=" << sizeof(Backend::value_type)
173 << "\n";
174 auto& prof = Alina::Profiler::globalProfiler();
175 Backend::params bprm;
176
177#if defined(SOLVER_BACKEND_CUDA)
178 cusparseCreate(&bprm.cusparse_handle);
179 {
180 int dev;
181 cudaGetDevice(&dev);
182
183 cudaDeviceProp prop;
184 cudaGetDeviceProperties(&prop, dev);
185 std::cout << prop.name << std::endl
186 << std::endl;
187 }
188#endif
189
191
193
194 if (reorder) {
195 prof.tic("reorder");
196 Alina::adapter::reorder<> perm(std::tie(rows, ptr, col, val));
197 prof.toc("reorder");
198
199 prof.tic("setup");
200 Solver solve(perm(std::tie(rows, ptr, col, val)), prm, bprm);
201 prof.toc("setup");
202
203 std::cout << solve << std::endl;
204
205 std::vector<double> tmp(rows);
206
207 perm.forward(rhs, tmp);
208 auto f_b = Backend::copy_vector(tmp, bprm);
209
210 perm.forward(x, tmp);
211 auto x_b = Backend::copy_vector(tmp, bprm);
212
213 prof.tic("solve");
214 info = solve(*f_b, *x_b);
215 prof.toc("solve");
216
217#if defined(SOLVER_BACKEND_CUDA)
218 thrust::copy(x_b->begin(), x_b->end(), tmp.begin());
219#else
220 std::copy(&(*x_b)[0], &(*x_b)[0] + rows, &tmp[0]);
221#endif
222
223 perm.inverse(tmp, x);
224 }
225 else {
226 prof.tic("setup");
227 Solver solve(std::tie(rows, ptr, col, val), prm, bprm);
228 prof.toc("setup");
229
230 std::cout << solve << std::endl;
231
232 auto f_b = Backend::copy_vector(rhs, bprm);
233 auto x_b = Backend::copy_vector(x, bprm);
234
235 prof.tic("solve");
236 info = solve(*f_b, *x_b);
237 prof.toc("solve");
238
239#if defined(SOLVER_BACKEND_CUDA)
240 thrust::copy(x_b->begin(), x_b->end(), x.begin());
241#else
242 std::copy(&(*x_b)[0], &(*x_b)[0] + rows, &x[0]);
243#endif
244 }
245
246 return info;
247}
248
249#define ARCCORE_ALINA_CALL_BLOCK_SOLVER(z, data, B) \
250 case B: \
251 return block_solve<B>(prm, rows, ptr, col, val, rhs, x, reorder);
252
253//---------------------------------------------------------------------------
255solve(const Alina::PropertyTree& prm,
256 size_t rows,
257 std::vector<ptrdiff_t> const& ptr,
258 std::vector<ptrdiff_t> const& col,
259 std::vector<double> const& val,
260 std::vector<double> const& rhs,
261 std::vector<double>& x,
262 int block_size,
263 bool reorder)
264{
265 switch (block_size) {
266 case 1:
267 return scalar_solve(prm, rows, ptr, col, val, rhs, x, reorder);
268#if defined(SOLVER_BACKEND_BUILTIN)
269 BOOST_PP_SEQ_FOR_EACH(ARCCORE_ALINA_CALL_BLOCK_SOLVER, ~, ARCCORE_ALINA_BLOCK_SIZES)
270#endif
271 default:
272 precondition(false, "Unsupported block size");
273 return {};
274 }
275}
276
277//---------------------------------------------------------------------------
278int main(int argc, char* argv[])
279{
280 auto& prof = Alina::Profiler::globalProfiler();
281 namespace po = boost::program_options;
282 namespace io = Alina::IO;
283
284 using std::string;
285 using std::vector;
286
287 po::options_description desc("Options");
288
289 desc.add_options()("help,h", "Show this help.")("prm-file,P",
290 po::value<string>(),
291 "Parameter file in json format. ")(
292 "prm,p",
293 po::value<vector<string>>()->multitoken(),
294 "Parameters specified as name=value pairs. "
295 "May be provided multiple times. Examples:\n"
296 " -p solver.tol=1e-3\n"
297 " -p precond.coarse_enough=300")("matrix,A",
298 po::value<string>(),
299 "System matrix in the MatrixMarket format. "
300 "When not specified, solves a Poisson problem in 3D unit cube. ")(
301 "rhs,f",
302 po::value<string>(),
303 "The RHS vector in the MatrixMarket format. "
304 "When omitted, a vector of ones is used by default. "
305 "Should only be provided together with a system matrix. ")(
306 "f0",
307 po::bool_switch()->default_value(false),
308 "Use zero RHS vector. Implies --random-initial and solver.ns_search=true")(
309 "f1",
310 po::bool_switch()->default_value(false),
311 "Set RHS = Ax where x = 1")(
312 "null,N",
313 po::value<string>(),
314 "The near null-space vectors in the MatrixMarket format. "
315 "Should be a dense matrix of size N*M, where N is the number of "
316 "unknowns, and M is the number of null-space vectors. "
317 "Should only be provided together with a system matrix. ")(
318 "coords,C",
319 po::value<string>(),
320 "Coordinate matrix where number of rows corresponds to the number of grid nodes "
321 "and the number of columns corresponds to the problem dimensionality (2 or 3). "
322 "Will be used to construct near null-space vectors as rigid body modes. "
323 "Should only be provided together with a system matrix. ")(
324 "binary,B",
325 po::bool_switch()->default_value(false),
326 "When specified, treat input files as binary instead of as MatrixMarket. "
327 "It is assumed the files were converted to binary format with mm2bin utility. ")(
328 "scale,s",
329 po::bool_switch()->default_value(false),
330 "Scale the matrix so that the diagonal is unit. ")(
331 "block-size,b",
332 po::value<int>()->default_value(1),
333 "The block size of the system matrix. "
334 "When specified, the system matrix is assumed to have block-wise structure. "
335 "This usually is the case for problems in elasticity, structural mechanics, "
336 "for coupled systems of PDE (such as Navier-Stokes equations), etc. ")(
337 "size,n",
338 po::value<int>()->default_value(32),
339 "The size of the Poisson problem to solve when no system matrix is given. "
340 "Specified as number of grid nodes along each dimension of a unit cube. "
341 "The resulting system will have n*n*n unknowns. ")(
342 "anisotropy,a",
343 po::value<double>()->default_value(1.0),
344 "The anisotropy value for the generated Poisson value. "
345 "Used to determine problem scaling along X, Y, and Z axes: "
346 "hy = hx * a, hz = hy * a.")(
347 "single-level,1",
348 po::bool_switch()->default_value(false),
349 "When specified, the AMG hierarchy is not constructed. "
350 "Instead, the problem is solved using a single-level smoother as preconditioner. ")(
351 "reorder,r",
352 po::bool_switch()->default_value(false),
353 "When specified, the matrix will be reordered to improve cache-locality")(
354 "initial,x",
355 po::value<double>()->default_value(0),
356 "Value to use as initial approximation. ")(
357 "random-initial",
358 po::bool_switch()->default_value(false),
359 "Use random initial approximation. ")(
360 "output,o",
361 po::value<string>(),
362 "Output file. Will be saved in the MatrixMarket format. "
363 "When omitted, the solution is not saved. ");
364
365 po::positional_options_description p;
366 p.add("prm", -1);
367
368 po::variables_map vm;
369 po::store(po::command_line_parser(argc, argv).options(desc).positional(p).run(), vm);
370 po::notify(vm);
371
372 if (vm.count("help")) {
373 std::cout << desc << std::endl;
374 return 0;
375 }
376
377 for (int i = 0; i < argc; ++i) {
378 if (i)
379 std::cout << " ";
380 std::cout << argv[i];
381 }
382 std::cout << std::endl;
383
385 if (vm.count("prm-file")) {
386 prm.read_json(vm["prm-file"].as<string>());
387 }
388
389 if (vm.count("prm")) {
390 for (const string& v : vm["prm"].as<vector<string>>()) {
391 prm.putKeyValue(v);
392 }
393 }
394
395 size_t rows, nv = 0;
396 vector<ptrdiff_t> ptr, col;
397 vector<double> val, rhs, null, x;
398
399 if (vm.count("matrix")) {
400 auto t = prof.scoped_tic("reading");
401
402 string Afile = vm["matrix"].as<string>();
403 bool binary = vm["binary"].as<bool>();
404
405 if (binary) {
406 io::read_crs(Afile, rows, ptr, col, val);
407 }
408 else {
409 size_t cols;
410 std::tie(rows, cols) = io::mm_reader(Afile)(ptr, col, val);
411 precondition(rows == cols, "Non-square system matrix");
412 }
413
414 if (vm.count("rhs")) {
415 string bfile = vm["rhs"].as<string>();
416
417 size_t n, m;
418
419 if (binary) {
420 io::read_dense(bfile, n, m, rhs);
421 }
422 else {
423 std::tie(n, m) = io::mm_reader(bfile)(rhs);
424 }
425
426 precondition(n == rows && m == 1, "The RHS vector has wrong size");
427 }
428 else if (vm["f1"].as<bool>()) {
429 rhs.resize(rows);
430 for (size_t i = 0; i < rows; ++i) {
431 double s = 0;
432 for (ptrdiff_t j = ptr[i], e = ptr[i + 1]; j < e; ++j)
433 s += val[j];
434 rhs[i] = s;
435 }
436 }
437 else {
438 rhs.resize(rows, vm["f0"].as<bool>() ? 0.0 : 1.0);
439 }
440
441 if (vm.count("null")) {
442 string nfile = vm["null"].as<string>();
443
444 size_t m;
445
446 if (binary) {
447 io::read_dense(nfile, m, nv, null);
448 }
449 else {
450 std::tie(m, nv) = io::mm_reader(nfile)(null);
451 }
452
453 precondition(m == rows, "Near null-space vectors have wrong size");
454 }
455 else if (vm.count("coords")) {
456 string cfile = vm["coords"].as<string>();
457 std::vector<double> coo;
458
459 size_t m, ndim;
460
461 if (binary) {
462 io::read_dense(cfile, m, ndim, coo);
463 }
464 else {
465 std::tie(m, ndim) = io::mm_reader(cfile)(coo);
466 }
467
468 precondition(m * ndim == rows && (ndim == 2 || ndim == 3), "Coordinate matrix has wrong size");
469
470 nv = Alina::rigid_body_modes(ndim, coo, null);
471 }
472
473 if (nv) {
474 prm.put("precond.coarsening.nullspace.cols", nv);
475 prm.put("precond.coarsening.nullspace.rows", rows);
476 prm.put("precond.coarsening.nullspace.B", &null[0]);
477 }
478 }
479 else {
480 auto t = prof.scoped_tic("assembling");
481 rows = sample_problem(vm["size"].as<int>(), val, col, ptr, rhs, vm["anisotropy"].as<double>());
482 }
483
484 if (vm["scale"].as<bool>()) {
485 std::vector<double> dia(rows, 1.0);
486
487 for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(rows); ++i) {
488 double d = 1.0;
489 for (ptrdiff_t j = ptr[i], e = ptr[i + 1]; j < e; ++j) {
490 if (col[j] == i) {
491 d = 1 / sqrt(val[j]);
492 }
493 }
494 if (!std::isnan(d))
495 dia[i] = d;
496 }
497
498 for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(rows); ++i) {
499 rhs[i] *= dia[i];
500 for (ptrdiff_t j = ptr[i], e = ptr[i + 1]; j < e; ++j) {
501 val[j] *= dia[i] * dia[col[j]];
502 }
503 }
504 }
505
506 x.resize(rows, vm["initial"].as<double>());
507 if (vm["random-initial"].as<bool>() || vm["f0"].as<bool>()) {
508 std::mt19937 rng;
509 std::uniform_real_distribution<double> rnd(-1, 1);
510 for (auto& v : x)
511 v = rnd(rng);
512 }
513
514 if (vm["f0"].as<bool>()) {
515 prm.put("solver.ns_search", true);
516 }
517
518 int block_size = vm["block-size"].as<int>();
519 std::cout << "BlockSize= " << block_size << "\n";
520
521 if (vm["single-level"].as<bool>())
522 prm.put("precond.class", "relaxation");
523
524 String do_hypre_str = Platform::getEnvironmentVariable("ALINA_USE_HYPRE");
525 bool do_hypre = false;
526 if (auto v = Convert::Type<Int32>::tryParseFromEnvironment("ALINA_USE_HYPRE", true))
527 do_hypre = v.value();
528
529 Alina::SolverResult solver_result;
530#ifndef SOLVER_BACKEND_BUILTIN
531 do_hypre = false;
532#endif
533 if (do_hypre) {
534#ifdef SOLVER_BACKEND_BUILTIN
535 _doHypreSolver(rows, ptr, col, val, rhs, x, argc, argv);
536#endif
537 }
538 else {
539 solver_result = solve(prm, rows, ptr, col, val, rhs, x, block_size, vm["reorder"].as<bool>());
540
541 if (vm.count("output")) {
542 auto t = prof.scoped_tic("write");
543 Alina::IO::mm_write(vm["output"].as<string>(), &x[0], x.size());
544 }
545 }
546 std::cout << "Iterations: " << solver_result.nbIteration() << std::endl
547 << "Error: " << solver_result.residual() << std::endl
548 << prof << std::endl;
549}
Convenience class that bundles together a preconditioner and an iterative solver.
Result of a solving.
Definition AlinaUtils.h:52
NUMA-aware vector container.
Definition NumaVector.h:42
Classe template pour convertir un type.
Vue d'un tableau d'éléments de type T.
Definition Span.h:801
Chaîne de caractères unicode.
ARCCORE_BASE_EXPORT String getEnvironmentVariable(const String &name)
Variable d'environnement du nom name.
__host__ __device__ double sqrt(double v)
Racine carrée de v.
Definition Math.h:135
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Alina::detail::empty_params params
Runtime-configurable wrappers around iterative solvers.