25#include <boost/program_options.hpp>
26#include <boost/range/iterator_range.hpp>
27#include <boost/scope_exit.hpp>
29#include "arccore/alina/BuiltinBackend.h"
30#include "arccore/alina/PreconditionerRuntime.h"
31#include "arccore/alina/Adapters.h"
32#include "arccore/alina/DistributedPreconditionedSolver.h"
33#include "arccore/alina/DistributedSolverRuntime.h"
34#include "arccore/alina/DistributedPreconditioner.h"
35#include "arccore/alina/Profiler.h"
36#include "arccore/alina/AlinaUtils.h"
39#include "arccore/alina/DistributedSolver.h"
42#include "arccore/alina/DistributedRelaxationRuntime.h"
44#include "DomainPartition.h"
47using namespace Arcane::Alina;
49using Alina::precondition;
54 const DomainPartition<2>& part;
55 const std::vector<ptrdiff_t>& dom;
57 renumbering(
const DomainPartition<2>& p,
58 const std::vector<ptrdiff_t>& d)
63 ptrdiff_t operator()(ptrdiff_t i, ptrdiff_t j)
const
65 boost::array<ptrdiff_t, 2> p = { { i, j } };
66 std::pair<int, ptrdiff_t> v = part.index(p);
67 return dom[v.first] + v.second;
73template <
template <
class>
class Precond,
class Matrix>
79 auto& prof = Alina::Profiler::globalProfiler();
85 const size_t n = Alina::backend::nbRow(A);
87 std::vector<double> rhs(n, 1), x(n, 0);
90 Solver solve(comm, A, prm);
94 auto t2 = prof.scoped_tic(
"solve");
101int main(
int argc,
char* argv[])
103 auto& prof = Alina::Profiler::globalProfiler();
104 namespace po = boost::program_options;
109 po::options_description desc(
"Options");
111 desc.add_options()(
"help,h",
"Show this help.")(
"prm-file,P",
113 "Parameter file in json format. ")(
115 po::value<vector<string>>()->multitoken(),
116 "Parameters specified as name=value pairs. "
117 "May be provided multiple times. Examples:\n"
118 " -p solver.tol=1e-3\n"
119 " -p precond.coarse_enough=300")(
121 po::value<int>()->default_value(1024),
122 "The size of the Poisson problem to solve. "
123 "Specified as number of grid nodes along each dimension of a unit square. "
124 "The resulting system will have n*n unknowns. ")(
126 po::bool_switch()->default_value(
false),
127 "When specified, the AMG hierarchy is not constructed. "
128 "Instead, the problem is solved using a single-level smoother as preconditioner. ")(
130 po::value<double>()->default_value(0),
131 "Value to use as initial approximation. ");
133 po::variables_map vm;
134 po::store(po::parse_command_line(argc, argv, desc), vm);
137 if (vm.count(
"help")) {
138 std::cout << desc << std::endl;
143 if (vm.count(
"prm-file")) {
144 prm.read_json(vm[
"prm-file"].as<string>());
147 if (vm.count(
"prm")) {
153 MPI_Init(&argc, &argv);
154 BOOST_SCOPE_EXIT(
void)
163 std::cout <<
"World size: " << world.size << std::endl;
165 const ptrdiff_t n = vm[
"size"].as<
int>();
166 const double h2i = (n - 1) * (n - 1);
168 boost::array<ptrdiff_t, 2> lo = { { 0, 0 } };
169 boost::array<ptrdiff_t, 2> hi = { { n - 1, n - 1 } };
171 prof.tic(
"partition");
172 DomainPartition<2> part(lo, hi, world.size);
173 ptrdiff_t chunk = part.size(world.rank);
175 std::vector<ptrdiff_t> domain(world.size + 1);
176 MPI_Allgather(&chunk, 1, Alina::mpi_datatype<ptrdiff_t>(),
177 &domain[1], 1, Alina::mpi_datatype<ptrdiff_t>(), world);
178 std::partial_sum(domain.begin(), domain.end(), domain.begin());
180 lo = part.domain(world.rank).min_corner();
181 hi = part.domain(world.rank).max_corner();
182 prof.toc(
"partition");
186 prof.tic(
"assemble");
187 std::vector<ptrdiff_t> ptr;
188 std::vector<ptrdiff_t> col;
189 std::vector<double> val;
190 std::vector<double> rhs;
192 ptr.reserve(chunk + 1);
193 col.reserve(chunk * 5);
194 val.reserve(chunk * 5);
198 for (ptrdiff_t j = lo[1]; j <= hi[1]; ++j) {
199 for (ptrdiff_t i = lo[0]; i <= hi[0]; ++i) {
201 col.push_back(renum(i, j - 1));
206 col.push_back(renum(i - 1, j));
210 col.push_back(renum(i, j));
211 val.push_back(4 * h2i);
214 col.push_back(renum(i + 1, j));
219 col.push_back(renum(i, j + 1));
223 ptr.push_back(col.size());
226 prof.toc(
"assemble");
228 bool single_level = vm[
"single-level"].as<
bool>();
231 prm.put(
"precond.class",
"relaxation");
233 Alina::SolverResult r = solve<Alina::PreconditionerRuntime>(world, prm, std::tie(chunk, ptr, col, val));
235 if (world.rank == 0) {
236 std::cout <<
"Iterations: " << r.nbIteration() << std::endl
237 <<
"Error: " << r.residual() << std::endl
239 << prof << std::endl;
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 -*-
Convenience wrapper around MPI_Comm.