Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
DistributedComplex.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 <complex>
23
24#include <boost/program_options.hpp>
25
26#include "arccore/alina/BuiltinBackend.h"
27#include "arccore/alina/ValueTypeComplex.h"
28#include "arccore/alina/Adapters.h"
29
30#include "arccore/alina/DistributedPreconditionedSolver.h"
31#include "arccore/alina/DistributedPreconditioner.h"
32#include "arccore/alina/DistributedSolverRuntime.h"
33
34#include "arccore/alina/IO.h"
35#include "arccore/alina/Profiler.h"
36
37using namespace Arcane;
38using namespace Arcane::Alina;
39
40namespace math = Alina::math;
41
42//---------------------------------------------------------------------------
43ptrdiff_t
44assemble_poisson3d(Alina::mpi_communicator comm,
45 ptrdiff_t n, int block_size,
46 std::vector<ptrdiff_t>& ptr,
47 std::vector<ptrdiff_t>& col,
48 std::vector<std::complex<double>>& val,
49 std::vector<std::complex<double>>& rhs)
50{
51 ptrdiff_t n3 = n * n * n;
52
53 ptrdiff_t chunk = (n3 + comm.size - 1) / comm.size;
54 if (chunk % block_size != 0) {
55 chunk += block_size - chunk % block_size;
56 }
57 ptrdiff_t row_beg = std::min(n3, chunk * comm.rank);
58 ptrdiff_t row_end = std::min(n3, row_beg + chunk);
59 chunk = row_end - row_beg;
60
61 ptr.clear();
62 ptr.reserve(chunk + 1);
63 col.clear();
64 col.reserve(chunk * 7);
65 val.clear();
66 val.reserve(chunk * 7);
67
68 rhs.resize(chunk);
69 std::fill(rhs.begin(), rhs.end(), 1.0);
70
71 const double h2i = (n - 1) * (n - 1);
72 ptr.push_back(0);
73
74 for (ptrdiff_t idx = row_beg; idx < row_end; ++idx) {
75 ptrdiff_t k = idx / (n * n);
76 ptrdiff_t j = (idx / n) % n;
77 ptrdiff_t i = idx % n;
78
79 if (k > 0) {
80 col.push_back(idx - n * n);
81 val.push_back(-h2i);
82 }
83
84 if (j > 0) {
85 col.push_back(idx - n);
86 val.push_back(-h2i);
87 }
88
89 if (i > 0) {
90 col.push_back(idx - 1);
91 val.push_back(-h2i);
92 }
93
94 col.push_back(idx);
95 val.push_back(6 * h2i);
96
97 if (i + 1 < n) {
98 col.push_back(idx + 1);
99 val.push_back(-h2i);
100 }
101
102 if (j + 1 < n) {
103 col.push_back(idx + n);
104 val.push_back(-h2i);
105 }
106
107 if (k + 1 < n) {
108 col.push_back(idx + n * n);
109 val.push_back(-h2i);
110 }
111
112 ptr.push_back(col.size());
113 }
114
115 return chunk;
116}
117
118//---------------------------------------------------------------------------
119void solve_scalar(Alina::mpi_communicator comm,
120 ptrdiff_t chunk,
121 const std::vector<ptrdiff_t>& ptr,
122 const std::vector<ptrdiff_t>& col,
123 const std::vector<std::complex<double>>& val,
124 const Alina::PropertyTree& prm,
125 const std::vector<std::complex<double>>& rhs)
126{
127 auto& prof = Alina::Profiler::globalProfiler();
129
132
133 prof.tic("setup");
134 Solver solve(comm, std::tie(chunk, ptr, col, val), prm);
135 prof.toc("setup");
136
137 if (comm.rank == 0) {
138 std::cout << solve << std::endl;
139 }
140
141 std::vector<std::complex<double>> x(chunk);
142
143 prof.tic("solve");
144 Alina::SolverResult r = solve(rhs, x);
145 prof.toc("solve");
146
147 if (comm.rank == 0) {
148 std::cout << "Iterations: " << r.nbIteration() << std::endl
149 << "Error: " << r.residual() << std::endl
150 << prof << std::endl;
151 }
152}
153
154//---------------------------------------------------------------------------
155int main(int argc, char* argv[])
156{
157 auto& prof = Alina::Profiler::globalProfiler();
158 Alina::mpi_init_thread mpi(&argc, &argv);
159 Alina::mpi_communicator comm(MPI_COMM_WORLD);
160
161 if (comm.rank == 0)
162 std::cout << "World size: " << comm.size << std::endl;
163
164 // Read configuration from command line
165 namespace po = boost::program_options;
166 po::options_description desc("Options");
167
168 desc.add_options()("help,h", "show help")(
169 "size,n",
170 po::value<ptrdiff_t>()->default_value(128),
171 "domain size")("prm-file,P",
172 po::value<std::string>(),
173 "Parameter file in json format. ")(
174 "prm,p",
175 po::value<std::vector<std::string>>()->multitoken(),
176 "Parameters specified as name=value pairs. "
177 "May be provided multiple times. Examples:\n"
178 " -p solver.tol=1e-3\n"
179 " -p precond.coarse_enough=300");
180
181 po::positional_options_description p;
182 p.add("prm", -1);
183
184 po::variables_map vm;
185 po::store(po::command_line_parser(argc, argv).options(desc).positional(p).run(), vm);
186 po::notify(vm);
187
188 if (vm.count("help")) {
189 if (comm.rank == 0)
190 std::cout << desc << std::endl;
191 return 0;
192 }
193
195 if (vm.count("prm-file")) {
196 prm.read_json(vm["prm-file"].as<std::string>());
197 }
198
199 if (vm.count("prm")) {
200 for (const std::string& v : vm["prm"].as<std::vector<std::string>>()) {
201 prm.putKeyValue(v);
202 }
203 }
204
205 ptrdiff_t n;
206 std::vector<ptrdiff_t> ptr;
207 std::vector<ptrdiff_t> col;
208 std::vector<std::complex<double>> val;
209 std::vector<std::complex<double>> rhs;
210
211 prof.tic("assemble");
212 n = assemble_poisson3d(comm, vm["size"].as<ptrdiff_t>(), 1, ptr, col, val, rhs);
213 prof.toc("assemble");
214
215 solve_scalar(comm, n, ptr, col, val, prm, rhs);
216}
Iterative solver wrapper for distributed linear systems.
Result of a solving.
Definition AlinaUtils.h:52
Espace de nom pour les fonctions mathématiques.
Definition MathUtils.h:41
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Convenience wrapper around MPI_Comm.
Convenience wrapper around MPI_Init_threads/MPI_Finalize.