20#pragma GCC diagnostic ignored "-Wdeprecated-copy"
21#pragma GCC diagnostic ignored "-Wint-in-bool-context"
23#include <boost/program_options.hpp>
24#include <boost/range/iterator_range.hpp>
25#include <boost/scope_exit.hpp>
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"
42using namespace Arcane::Alina;
44using Alina::precondition;
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)
56 ptrdiff_t n = A_mm.rows();
58 ptrdiff_t chunk = (n + comm.size - 1) / comm.size;
59 if (chunk % block_size != 0) {
60 chunk += block_size - chunk % block_size;
63 ptrdiff_t row_beg = std::min(n, chunk * comm.rank);
64 ptrdiff_t row_end = std::min(n, row_beg + chunk);
66 chunk = row_end - row_beg;
68 A_mm(ptr, col, val, row_beg, row_end);
70 if (rhs_file.empty()) {
72 std::fill(rhs.begin(), rhs.end(), 1.0);
76 rhs_mm(rhs, row_beg, row_end);
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)
91 ptrdiff_t n = Alina::IO::crs_size<ptrdiff_t>(A_file);
93 ptrdiff_t chunk = (n + comm.size - 1) / comm.size;
94 if (chunk % block_size != 0) {
95 chunk += block_size - chunk % block_size;
98 ptrdiff_t row_beg = std::min(n, chunk * comm.rank);
99 ptrdiff_t row_end = std::min(n, row_beg + chunk);
101 chunk = row_end - row_beg;
103 Alina::IO::read_crs(A_file, n, ptr, col, val, row_beg, row_end);
105 if (rhs_file.empty()) {
107 std::fill(rhs.begin(), rhs.end(), 1.0);
110 ptrdiff_t rows, cols;
111 Alina::IO::read_dense(rhs_file, rows, cols, rhs, row_beg, row_end);
118template <
class Backend,
class Matrix>
119std::shared_ptr<Alina::DistributedMatrix<Backend>>
122 Alina::eMatrixPartitionerType ptype,
int block_size = 1)
124 auto& prof = Alina::Profiler::globalProfiler();
127 auto A = std::make_shared<DMatrix>(comm, Astrip);
129 if (comm.size == 1 || ptype == Alina::eMatrixPartitionerType::merge)
132 prof.tic(
"partition");
134 prm.put(
"type", ptype);
137 auto I = part(*A, block_size);
138 auto J = transpose(*I);
139 A = product(*J, *product(*A, *I));
141 std::vector<double> new_rhs(J->loc_rows());
143 J->move_to_backend(bprm);
145 Alina::backend::spmv(1, *J, rhs, 0, new_rhs);
147 prof.toc(
"partition");
154int main(
int argc,
char* argv[])
156 auto& prof = Alina::Profiler::globalProfiler();
158 MPI_Init_thread(&argc, &argv, MPI_THREAD_MULTIPLE, &provided);
159 BOOST_SCOPE_EXIT(
void)
168 std::cout <<
"World size: " << comm.size << std::endl;
171 namespace po = boost::program_options;
172 po::options_description desc(
"Options");
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. ")(
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. ")(
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. ")(
188 po::value<int>()->default_value(1),
189 "The block size of the system matrix. ")(
191 po::value<Alina::eMatrixPartitionerType>()->default_value(
192#
if defined(ARCCORE_ALINA_HAVE_PARMETIS)
193 Alina::eMatrixPartitionerType::parmetis
196 "Repartition the system matrix")(
"prm-file,P",
197 po::value<std::string>(),
198 "Parameter file in json format. ")(
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");
206 po::variables_map vm;
207 po::store(po::parse_command_line(argc, argv, desc), vm);
210 if (vm.count(
"help")) {
212 std::cout << desc << std::endl;
217 if (vm.count(
"prm-file")) {
218 prm.read_json(vm[
"prm-file"].as<std::string>());
221 if (vm.count(
"prm")) {
222 for (
const std::string& v : vm[
"prm"].as<std::vector<std::string>>()) {
228 std::vector<ptrdiff_t> ptr;
229 std::vector<ptrdiff_t> col;
230 std::vector<double> val;
231 std::vector<double> rhs;
233 int block_size = vm[
"block-size"].as<
int>();
234 prm.put(
"precond.block_size", block_size);
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);
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);
253 auto A = partition<Backend>(comm,
255 vm[
"partitioner"].as<Alina::eMatrixPartitionerType>(),
273 Solver solve(comm, A, prm);
277 std::cout << solve << std::endl;
279 std::vector<double> x(rhs.size(), 0.0);
285 if (comm.rank == 0) {
286 std::cout <<
"Iterations: " << r.nbIteration() << std::endl
287 <<
"Error: " << r.residual() << std::endl
288 << prof << std::endl;
Algebraic multigrid method.
Distributed CPR preconditioner.
Runtime wrapper for distributed direct solvers.
Distributed Matrix using message passing.
Iterative solver wrapper for distributed linear systems.
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.