Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
DistributedCPR.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// To remove warnings about deprecated Eigen usage.
20#pragma GCC diagnostic ignored "-Wdeprecated-copy"
21#pragma GCC diagnostic ignored "-Wint-in-bool-context"
22
23#include <boost/program_options.hpp>
24#include <boost/range/iterator_range.hpp>
25#include <boost/scope_exit.hpp>
26
27#include "arccore/alina/IO.h"
28#include "arccore/alina/Adapters.h"
29#include "arccore/alina/BuiltinBackend.h"
30#include "arccore/alina/DistributedPreconditionedSolver.h"
31#include "arccore/alina/DistributedCPRPreconditioner.h"
32#include "arccore/alina/DistributedAMG.h"
33#include "arccore/alina/DistributedCoarseningRuntime.h"
34#include "arccore/alina/DistributedRelaxationRuntime.h"
35#include "arccore/alina/DistributedSolverRuntime.h"
36#include "arccore/alina/DistributedDirectSolverRuntime.h"
37#include "arccore/alina/MatrixPartitionerRuntime.h"
38#include "arccore/alina/Profiler.h"
39#include "arccore/alina/AlinaUtils.h"
40
41using namespace Arcane;
42using namespace Arcane::Alina;
43
44using Alina::precondition;
45
46//---------------------------------------------------------------------------
47ptrdiff_t
48read_matrix_market(Alina::mpi_communicator comm,
49 const std::string& A_file, const std::string& rhs_file, int block_size,
50 std::vector<ptrdiff_t>& ptr,
51 std::vector<ptrdiff_t>& col,
52 std::vector<double>& val,
53 std::vector<double>& rhs)
54{
55 Alina::IO::mm_reader A_mm(A_file);
56 ptrdiff_t n = A_mm.rows();
57
58 ptrdiff_t chunk = (n + comm.size - 1) / comm.size;
59 if (chunk % block_size != 0) {
60 chunk += block_size - chunk % block_size;
61 }
62
63 ptrdiff_t row_beg = std::min(n, chunk * comm.rank);
64 ptrdiff_t row_end = std::min(n, row_beg + chunk);
65
66 chunk = row_end - row_beg;
67
68 A_mm(ptr, col, val, row_beg, row_end);
69
70 if (rhs_file.empty()) {
71 rhs.resize(chunk);
72 std::fill(rhs.begin(), rhs.end(), 1.0);
73 }
74 else {
75 Alina::IO::mm_reader rhs_mm(rhs_file);
76 rhs_mm(rhs, row_beg, row_end);
77 }
78
79 return chunk;
80}
81
82//---------------------------------------------------------------------------
83ptrdiff_t
84read_binary(Alina::mpi_communicator comm,
85 const std::string& A_file, const std::string& rhs_file, int block_size,
86 std::vector<ptrdiff_t>& ptr,
87 std::vector<ptrdiff_t>& col,
88 std::vector<double>& val,
89 std::vector<double>& rhs)
90{
91 ptrdiff_t n = Alina::IO::crs_size<ptrdiff_t>(A_file);
92
93 ptrdiff_t chunk = (n + comm.size - 1) / comm.size;
94 if (chunk % block_size != 0) {
95 chunk += block_size - chunk % block_size;
96 }
97
98 ptrdiff_t row_beg = std::min(n, chunk * comm.rank);
99 ptrdiff_t row_end = std::min(n, row_beg + chunk);
100
101 chunk = row_end - row_beg;
102
103 Alina::IO::read_crs(A_file, n, ptr, col, val, row_beg, row_end);
104
105 if (rhs_file.empty()) {
106 rhs.resize(chunk);
107 std::fill(rhs.begin(), rhs.end(), 1.0);
108 }
109 else {
110 ptrdiff_t rows, cols;
111 Alina::IO::read_dense(rhs_file, rows, cols, rhs, row_beg, row_end);
112 }
113
114 return chunk;
115}
116
117//---------------------------------------------------------------------------
118template <class Backend, class Matrix>
119std::shared_ptr<Alina::DistributedMatrix<Backend>>
120partition(Alina::mpi_communicator comm, const Matrix& Astrip,
121 std::vector<double>& rhs, const typename Backend::params& bprm,
122 Alina::eMatrixPartitionerType ptype, int block_size = 1)
123{
124 auto& prof = Alina::Profiler::globalProfiler();
125 typedef Alina::DistributedMatrix<Backend> DMatrix;
126
127 auto A = std::make_shared<DMatrix>(comm, Astrip);
128
129 if (comm.size == 1 || ptype == Alina::eMatrixPartitionerType::merge)
130 return A;
131
132 prof.tic("partition");
134 prm.put("type", ptype);
136
137 auto I = part(*A, block_size);
138 auto J = transpose(*I);
139 A = product(*J, *product(*A, *I));
140
141 std::vector<double> new_rhs(J->loc_rows());
142
143 J->move_to_backend(bprm);
144
145 Alina::backend::spmv(1, *J, rhs, 0, new_rhs);
146 rhs.swap(new_rhs);
147 prof.toc("partition");
148
149 return A;
150}
151
152//---------------------------------------------------------------------------
153
154int main(int argc, char* argv[])
155{
156 auto& prof = Alina::Profiler::globalProfiler();
157 int provided;
158 MPI_Init_thread(&argc, &argv, MPI_THREAD_MULTIPLE, &provided);
159 BOOST_SCOPE_EXIT(void)
160 {
161 MPI_Finalize();
162 }
163 BOOST_SCOPE_EXIT_END
164
165 Alina::mpi_communicator comm(MPI_COMM_WORLD);
166
167 if (comm.rank == 0)
168 std::cout << "World size: " << comm.size << std::endl;
169
170 // Read configuration from command line
171 namespace po = boost::program_options;
172 po::options_description desc("Options");
173
174 desc.add_options()("help,h", "show help")("matrix,A",
175 po::value<std::string>(),
176 "System matrix in the MatrixMarket format. "
177 "When not specified, a Poisson problem in 3D unit cube is assembled. ")(
178 "rhs,f",
179 po::value<std::string>()->default_value(""),
180 "The RHS vector in the MatrixMarket format. "
181 "When omitted, a vector of ones is used by default. "
182 "Should only be provided together with a system matrix. ")(
183 "binary,B",
184 po::bool_switch()->default_value(false),
185 "When specified, treat input files as binary instead of as MatrixMarket. "
186 "It is assumed the files were converted to binary format with mm2bin utility. ")(
187 "block-size,b",
188 po::value<int>()->default_value(1),
189 "The block size of the system matrix. ")(
190 "partitioner,r",
191 po::value<Alina::eMatrixPartitionerType>()->default_value(
192#if defined(ARCCORE_ALINA_HAVE_PARMETIS)
193 Alina::eMatrixPartitionerType::parmetis
194#endif
195 ),
196 "Repartition the system matrix")("prm-file,P",
197 po::value<std::string>(),
198 "Parameter file in json format. ")(
199 "prm,p",
200 po::value<std::vector<std::string>>()->multitoken(),
201 "Parameters specified as name=value pairs. "
202 "May be provided multiple times. Examples:\n"
203 " -p solver.tol=1e-3\n"
204 " -p precond.coarse_enough=300");
205
206 po::variables_map vm;
207 po::store(po::parse_command_line(argc, argv, desc), vm);
208 po::notify(vm);
209
210 if (vm.count("help")) {
211 if (comm.rank == 0)
212 std::cout << desc << std::endl;
213 return 0;
214 }
215
217 if (vm.count("prm-file")) {
218 prm.read_json(vm["prm-file"].as<std::string>());
219 }
220
221 if (vm.count("prm")) {
222 for (const std::string& v : vm["prm"].as<std::vector<std::string>>()) {
223 prm.putKeyValue(v);
224 }
225 }
226
227 ptrdiff_t n;
228 std::vector<ptrdiff_t> ptr;
229 std::vector<ptrdiff_t> col;
230 std::vector<double> val;
231 std::vector<double> rhs;
232
233 int block_size = vm["block-size"].as<int>();
234 prm.put("precond.block_size", block_size);
235
236 prof.tic("read");
237 if (vm["binary"].as<bool>()) {
238 n = read_binary(comm,
239 vm["matrix"].as<std::string>(),
240 vm["rhs"].as<std::string>(),
241 block_size, ptr, col, val, rhs);
242 }
243 else {
244 n = read_matrix_market(comm,
245 vm["matrix"].as<std::string>(),
246 vm["rhs"].as<std::string>(),
247 block_size, ptr, col, val, rhs);
248 }
249 prof.toc("read");
250
251 typedef Alina::BuiltinBackend<double> Backend;
252
253 auto A = partition<Backend>(comm,
254 std::tie(n, ptr, col, val), rhs, Backend::params(),
255 vm["partitioner"].as<Alina::eMatrixPartitionerType>(),
256 block_size);
257
258 prof.tic("setup");
259
260 using AMG = DistributedAMG<Backend,
265
268 AMG,
271 Solver;
272
273 Solver solve(comm, A, prm);
274 prof.toc("setup");
275
276 if (comm.rank == 0)
277 std::cout << solve << std::endl;
278
279 std::vector<double> x(rhs.size(), 0.0);
280
281 prof.tic("solve");
282 Alina::SolverResult r = solve(rhs, x);
283 prof.toc("solve");
284
285 if (comm.rank == 0) {
286 std::cout << "Iterations: " << r.nbIteration() << std::endl
287 << "Error: " << r.residual() << std::endl
288 << prof << std::endl;
289 }
290}
Algebraic multigrid method.
Definition AMG.h:71
Runtime wrapper for distributed direct solvers.
Distributed Matrix using message passing.
Iterative solver wrapper for distributed linear systems.
Matrix market reader.
Definition IO.h:52
Result of a solving.
Definition AlinaUtils.h:52
Matrix class, to be used by user.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Use a relaxation as a distributed preconditioner.
Alina::detail::empty_params params
Distributed memory sparse approximate inverse relaxation scheme.
Runtime-configurable wrapper around matrix partitioner.
Convenience wrapper around MPI_Comm.