Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
DistributedRuntimeBP.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 <vector>
21#include <string>
22#include <utility>
23#include <numeric>
24
25#include <boost/program_options.hpp>
26#include <boost/range/iterator_range.hpp>
27#include <boost/scope_exit.hpp>
28
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"
37
38// Pour test compilation uniquement
39#include "arccore/alina/DistributedSolver.h"
40
41// Pour test compilation uniquement
42#include "arccore/alina/DistributedRelaxationRuntime.h"
43
44#include "DomainPartition.h"
45
46using namespace Arcane;
47using namespace Arcane::Alina;
48
49using Alina::precondition;
50
51//---------------------------------------------------------------------------
52struct renumbering
53{
54 const DomainPartition<2>& part;
55 const std::vector<ptrdiff_t>& dom;
56
57 renumbering(const DomainPartition<2>& p,
58 const std::vector<ptrdiff_t>& d)
59 : part(p)
60 , dom(d)
61 {}
62
63 ptrdiff_t operator()(ptrdiff_t i, ptrdiff_t j) const
64 {
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;
68 }
69};
70
71//---------------------------------------------------------------------------
72
73template <template <class> class Precond, class Matrix>
75solve(const Alina::mpi_communicator& comm,
76 const Alina::PropertyTree& prm,
77 const Matrix& A)
78{
79 auto& prof = Alina::Profiler::globalProfiler();
80
81 typedef Alina::BuiltinBackend<double> Backend;
82
84
85 const size_t n = Alina::backend::nbRow(A);
86
87 std::vector<double> rhs(n, 1), x(n, 0);
88
89 prof.tic("setup");
90 Solver solve(comm, A, prm);
91 prof.toc("setup");
92
93 {
94 auto t2 = prof.scoped_tic("solve");
95 return solve(rhs, x);
96 }
97}
98
99//---------------------------------------------------------------------------
100
101int main(int argc, char* argv[])
102{
103 auto& prof = Alina::Profiler::globalProfiler();
104 namespace po = boost::program_options;
105
106 using std::string;
107 using std::vector;
108
109 po::options_description desc("Options");
110
111 desc.add_options()("help,h", "Show this help.")("prm-file,P",
112 po::value<string>(),
113 "Parameter file in json format. ")(
114 "prm,p",
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")(
120 "size,n",
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. ")(
125 "single-level,1",
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. ")(
129 "initial,x",
130 po::value<double>()->default_value(0),
131 "Value to use as initial approximation. ");
132
133 po::variables_map vm;
134 po::store(po::parse_command_line(argc, argv, desc), vm);
135 po::notify(vm);
136
137 if (vm.count("help")) {
138 std::cout << desc << std::endl;
139 return 0;
140 }
141
143 if (vm.count("prm-file")) {
144 prm.read_json(vm["prm-file"].as<string>());
145 }
146
147 if (vm.count("prm")) {
148 for (const string& v : vm["prm"].as<vector<string>>()) {
149 prm.putKeyValue(v);
150 }
151 }
152
153 MPI_Init(&argc, &argv);
154 BOOST_SCOPE_EXIT(void)
155 {
156 MPI_Finalize();
157 }
158 BOOST_SCOPE_EXIT_END
159
160 Alina::mpi_communicator world(MPI_COMM_WORLD);
161
162 if (world.rank == 0)
163 std::cout << "World size: " << world.size << std::endl;
164
165 const ptrdiff_t n = vm["size"].as<int>();
166 const double h2i = (n - 1) * (n - 1);
167
168 boost::array<ptrdiff_t, 2> lo = { { 0, 0 } };
169 boost::array<ptrdiff_t, 2> hi = { { n - 1, n - 1 } };
170
171 prof.tic("partition");
172 DomainPartition<2> part(lo, hi, world.size);
173 ptrdiff_t chunk = part.size(world.rank);
174
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());
179
180 lo = part.domain(world.rank).min_corner();
181 hi = part.domain(world.rank).max_corner();
182 prof.toc("partition");
183
184 renumbering renum(part, domain);
185
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;
191
192 ptr.reserve(chunk + 1);
193 col.reserve(chunk * 5);
194 val.reserve(chunk * 5);
195
196 ptr.push_back(0);
197
198 for (ptrdiff_t j = lo[1]; j <= hi[1]; ++j) {
199 for (ptrdiff_t i = lo[0]; i <= hi[0]; ++i) {
200 if (j > 0) {
201 col.push_back(renum(i, j - 1));
202 val.push_back(-h2i);
203 }
204
205 if (i > 0) {
206 col.push_back(renum(i - 1, j));
207 val.push_back(-h2i);
208 }
209
210 col.push_back(renum(i, j));
211 val.push_back(4 * h2i);
212
213 if (i + 1 < n) {
214 col.push_back(renum(i + 1, j));
215 val.push_back(-h2i);
216 }
217
218 if (j + 1 < n) {
219 col.push_back(renum(i, j + 1));
220 val.push_back(-h2i);
221 }
222
223 ptr.push_back(col.size());
224 }
225 }
226 prof.toc("assemble");
227
228 bool single_level = vm["single-level"].as<bool>();
229
230 if (single_level)
231 prm.put("precond.class", "relaxation");
232
233 Alina::SolverResult r = solve<Alina::PreconditionerRuntime>(world, prm, std::tie(chunk, ptr, col, val));
234
235 if (world.rank == 0) {
236 std::cout << "Iterations: " << r.nbIteration() << std::endl
237 << "Error: " << r.residual() << std::endl
238 << std::endl
239 << prof << std::endl;
240 }
241}
Iterative solver wrapper for distributed linear systems.
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 -*-
Convenience wrapper around MPI_Comm.