23#include <boost/program_options.hpp>
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"
33#include "arccore/alina/IO.h"
34#include "arccore/alina/Profiler.h"
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"
43#include "AlinaSamplesCommon.h"
47namespace math = Alina::math;
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)
58 ptrdiff_t n3 = n * n * n;
60 ptrdiff_t chunk = (n3 + comm.size - 1) / comm.size;
61 if (chunk % block_size != 0) {
62 chunk += block_size - chunk % block_size;
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;
69 ptr.reserve(chunk + 1);
71 col.reserve(chunk * 7);
73 val.reserve(chunk * 7);
76 std::fill(rhs.begin(), rhs.end(), 1.0);
78 const double h2i = (n - 1) * (n - 1);
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;
87 col.push_back(idx - n * n);
92 col.push_back(idx - n);
97 col.push_back(idx - 1);
102 val.push_back(6 * h2i);
105 col.push_back(idx + 1);
110 col.push_back(idx + n);
115 col.push_back(idx + n * n);
119 ptr.push_back(col.size());
126template <
class Backend,
class Matrix>
127std::shared_ptr<Alina::DistributedMatrix<Backend>>
130 Alina::eMatrixPartitionerType ptype,
int block_size = 1)
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;
137 auto A = std::make_shared<DMatrix>(comm, Astrip);
139 if (comm.size == 1 || ptype == Alina::eMatrixPartitionerType::merge)
142 prof.tic(
"partition");
144 prm.put(
"type", ptype);
147 auto I = part(*A, block_size);
149 A = product(*J, *product(*A, *I));
153 J->move_to_backend(bprm);
155 Alina::backend::spmv(1, *J, rhs, 0, new_rhs);
157 prof.toc(
"partition");
165 const std::vector<ptrdiff_t>& ptr,
166 const std::vector<ptrdiff_t>& col,
167 const std::vector<double>& val,
169 const std::vector<double>& f,
170 Alina::eMatrixPartitionerType ptype)
172 auto& prof = Alina::Profiler::globalProfiler();
175 using BackendValueType = double;
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)
202 auto get_distributed_matrix = [&]() {
203 auto t = prof.scoped_tic(
"distributed matrix");
204 std::shared_ptr<DMatrix> A;
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();
213 A = std::make_shared<DMatrix>(comm, std::tie(chunk, ptr, col, val));
219 std::shared_ptr<DMatrix> A;
220 std::shared_ptr<Solver> solve;
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";
231 if (comm.rank == 0) {
232 std::cout <<
"SolverInfo:\n";
233 std::cout << *solve << std::endl;
236 if (prm.get(
"precond.allow_rebuild",
false)) {
237 if (comm.rank == 0) {
238 std::cout <<
"Rebuilding the preconditioner..." << std::endl;
242 auto t = prof.scoped_tic(
"rebuild");
243 A = get_distributed_matrix();
244 solve->precond().rebuild(A);
247 if (comm.rank == 0) {
248 std::cout << *solve << std::endl;
258 if (comm.rank == 0) {
259 std::cout <<
"Iterations: " << r.nbIteration() << std::endl
260 <<
"Error: " << r.residual() << std::endl
261 << prof << std::endl;
269 auto& prof = Alina::Profiler::globalProfiler();
274 tm->
info() <<
"World size: " << comm.size;
277 namespace po = boost::program_options;
278 po::options_description desc(
"Options");
280 auto default_partitioner_type = Alina::eMatrixPartitionerType::merge;
281#if defined(ARCCORE_ALINA_HAVE_PARMETIS)
282 default_partitioner_type = Alina::eMatrixPartitionerType::parmetis;
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),
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. ");
308 po::positional_options_description p;
311 po::variables_map vm;
312 po::store(po::command_line_parser(argc, argv).options(desc).positional(p).run(), vm);
315 if (vm.count(
"help")) {
317 std::cout << desc << std::endl;
322 if (vm.count(
"prm-file")) {
323 prm.read_json(vm[
"prm-file"].as<std::string>());
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;
334 std::vector<ptrdiff_t> ptr;
335 std::vector<ptrdiff_t> col;
336 std::vector<double> val;
337 std::vector<double> rhs;
339 Alina::eMatrixPartitionerType ptype = vm[
"partitioner"].as<Alina::eMatrixPartitionerType>();
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");
347 if (vm[
"test-rebuild"].as<bool>()) {
348 prm.put(
"precond.allow_rebuild",
true);
351 solve_scalar(comm, n, ptr, col, val, prm, rhs, ptype);
358int main(
int argc,
char* argv[])
360 return Arcane::Alina::SampleMainContext::execMain(main2, argc, argv);
Runtime wrapper for distributed direct solvers.
Distributed Matrix using message passing.
Iterative solver wrapper for distributed linear systems.
NUMA-aware vector container.
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.
__host__ __device__ Real3x3 transpose(const Real3x3 &t)
Transpose la matrice.
-*- 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.