Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
CPRDynamicRowSum.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 <string>
21
22#include <boost/program_options.hpp>
23#include <boost/preprocessor/seq/for_each.hpp>
24
25#if defined(SOLVER_BACKEND_CUDA)
26# include "arccore/alina/CudaBackend.h"
27# include "arccore/alina/relaxation_cusparse_ilu0.h"
29#else
30# ifndef SOLVER_BACKEND_BUILTIN
31# define SOLVER_BACKEND_BUILTIN
32# endif
33#include "arccore/alina/BuiltinBackend.h"
35#endif
36
37#if defined(SOLVER_BACKEND_BUILTIN)
38#include "arccore/alina/StaticMatrix.h"
39#include "arccore/alina/Adapters.h"
40#endif
41
42#include "arccore/alina/PreconditionedSolver.h"
43#include "arccore/alina/AMG.h"
44#include "arccore/alina/SolverRuntime.h"
45#include "arccore/alina/CoarseningRuntime.h"
46#include "arccore/alina/RelaxationRuntime.h"
47#include "arccore/alina/Relaxation.h"
48#include "arccore/alina/CPRDynamicRowSumPreconditioner.h"
49#include "arccore/alina/Adapters.h"
50#include "arccore/alina/IO.h"
51#include "arccore/alina/Profiler.h"
52
53using namespace Arcane;
54
55using Alina::precondition;
56
57//---------------------------------------------------------------------------
58template <class Matrix>
59void solve_cpr(const Matrix& K, const std::vector<double>& rhs, Alina::PropertyTree& prm)
60{
61 auto& prof = Alina::Profiler::globalProfiler();
62 Backend::params bprm;
63
64#if defined(SOLVER_BACKEND_VEXCL)
65 vex::Context ctx(vex::Filter::Env);
66 std::cout << ctx << std::endl;
67 bprm.q = ctx;
68#elif defined(SOLVER_BACKEND_VIENNACL)
69 std::cout
70 << viennacl::ocl::current_device().name()
71 << " (" << viennacl::ocl::current_device().vendor() << ")\n\n";
72#elif defined(SOLVER_BACKEND_CUDA)
73 cusparseCreate(&bprm.cusparse_handle);
74 {
75 int dev;
76 cudaGetDevice(&dev);
77
78 cudaDeviceProp prop;
79 cudaGetDeviceProperties(&prop, dev);
80 std::cout << prop.name << std::endl
81 << std::endl;
82 }
83#endif
84
85 auto t1 = prof.scoped_tic("CPR");
86
88 PPrecond;
89
91
92 prof.tic("setup");
96 solve(K, prm, bprm);
97 prof.toc("setup");
98
99 std::cout << solve.precond() << std::endl;
100
101 auto f = Backend::copy_vector(rhs, bprm);
102 auto x = Backend::create_vector(rhs.size(), bprm);
103 Alina::backend::clear(*x);
104
105 prof.tic("solve");
106 Alina::SolverResult r = solve(*f, *x);
107 prof.toc("solve");
108
109 std::cout << "Iterations: " << r.nbIteration() << std::endl
110 << "Error: " << r.residual() << std::endl;
111}
112
113#if defined(SOLVER_BACKEND_BUILTIN)
114//---------------------------------------------------------------------------
115template <int B, class Matrix>
116void solve_block_cpr(const Matrix& K, const std::vector<double>& rhs, Alina::PropertyTree& prm)
117{
118 auto& prof = Alina::Profiler::globalProfiler();
119 auto t1 = prof.scoped_tic("CPR");
120
121 typedef Alina::StaticMatrix<double, B, B> val_type;
122 typedef Alina::StaticMatrix<double, B, 1> rhs_type;
123
124 typedef Alina::BuiltinBackend<val_type> SBackend;
125 typedef Alina::BuiltinBackend<double> PBackend;
126
129
130 typename SBackend::params bprm;
131
132 prof.tic("setup");
136 solve(Alina::adapter::block_matrix<val_type>(K), prm, bprm);
137 prof.toc("setup");
138
139 std::cout << solve.precond() << std::endl;
140
141 size_t n = Alina::backend::nbRow(K) / B;
142 auto rhs_ptr = reinterpret_cast<const rhs_type*>(rhs.data());
143
144 SmallSpan<const rhs_type> f(rhs_ptr, n);
145
146 auto x = SBackend::create_vector(n, bprm);
147 Alina::backend::clear(*x);
148
149 prof.tic("solve");
150 Alina::SolverResult r = solve(f, *x);
151 prof.toc("solve");
152
153 std::cout << "Iterations: " << r.nbIteration() << std::endl
154 << "Error: " << r.residual() << std::endl;
155}
156#endif
157
158//---------------------------------------------------------------------------
159int main(int argc, char* argv[])
160{
161 auto& prof = Alina::Profiler::globalProfiler();
162 using Alina::precondition;
163 using std::string;
164 using std::vector;
165
166 namespace po = boost::program_options;
167 namespace io = Alina::IO;
168
169 po::options_description desc("Options");
170
171 desc.add_options()("help,h", "show help")(
172 "binary,B",
173 po::bool_switch()->default_value(false),
174 "When specified, treat input files as binary instead of as MatrixMarket. "
175 "It is assumed the files were converted to binary format with mm2bin utility. ")(
176 "matrix,A",
177 po::value<string>()->required(),
178 "The system matrix in MatrixMarket format")(
179 "rhs,f",
180 po::value<string>(),
181 "The right-hand side in MatrixMarket format")(
182 "weights,w",
183 po::value<string>(),
184 "Equation weights in MatrixMarket format")(
185 "runtime-block-size,b",
186 po::value<int>(),
187 "The block size of the system matrix set at runtime")(
188 "static-block-size,c",
189 po::value<int>()->default_value(1),
190 "The block size of the system matrix set at compiletime")(
191 "params,P",
192 po::value<string>(),
193 "parameter file in json format")(
194 "prm,p",
195 po::value<vector<string>>()->multitoken(),
196 "Parameters specified as name=value pairs. "
197 "May be provided multiple times. Examples:\n"
198 " -p solver.tol=1e-3\n"
199 " -p precond.coarse_enough=300");
200
201 po::variables_map vm;
202 po::store(po::parse_command_line(argc, argv, desc), vm);
203
204 if (vm.count("help")) {
205 std::cout << desc << std::endl;
206 return 0;
207 }
208
209 po::notify(vm);
210
212 if (vm.count("params"))
213 prm.read_json(vm["params"].as<string>());
214
215 if (vm.count("prm")) {
216 for (const string& v : vm["prm"].as<vector<string>>()) {
217 prm.putKeyValue(v);
218 }
219 }
220
221 int cb = vm["static-block-size"].as<int>();
222
223 if (vm.count("runtime-block-size"))
224 prm.put("precond.block_size", vm["runtime-block-size"].as<int>());
225 else
226 prm.put("precond.block_size", cb);
227
228 size_t rows;
229 vector<ptrdiff_t> ptr, col;
230 vector<double> val, rhs, wgt;
231 std::vector<char> pm;
232
233 {
234 auto t = prof.scoped_tic("reading");
235
236 string Afile = vm["matrix"].as<string>();
237 bool binary = vm["binary"].as<bool>();
238
239 if (binary) {
240 io::read_crs(Afile, rows, ptr, col, val);
241 }
242 else {
243 size_t cols;
244 std::tie(rows, cols) = io::mm_reader(Afile)(ptr, col, val);
245 precondition(rows == cols, "Non-square system matrix");
246 }
247
248 if (vm.count("rhs")) {
249 string bfile = vm["rhs"].as<string>();
250
251 size_t n, m;
252
253 if (binary) {
254 io::read_dense(bfile, n, m, rhs);
255 }
256 else {
257 std::tie(n, m) = io::mm_reader(bfile)(rhs);
258 }
259
260 precondition(n == rows && m == 1, "The RHS vector has wrong size");
261 }
262 else {
263 rhs.resize(rows, 1.0);
264 }
265
266 if (vm.count("weights")) {
267 string wfile = vm["weights"].as<string>();
268
269 size_t n, m;
270
271 if (binary) {
272 io::read_dense(wfile, n, m, wgt);
273 }
274 else {
275 std::tie(n, m) = io::mm_reader(wfile)(wgt);
276 }
277
278 prm.put("precond.weights", &wgt[0]);
279 prm.put("precond.weights_size", wgt.size());
280 }
281 }
282
283#define CALL_BLOCK_SOLVER(z, data, B) \
284 case B: \
285 solve_block_cpr<B>(std::tie(rows, ptr, col, val), rhs, prm); \
286 break;
287
288 switch (cb) {
289 case 1:
290 solve_cpr(std::tie(rows, ptr, col, val), rhs, prm);
291 break;
292
293#if defined(SOLVER_BACKEND_BUILTIN) || defined(SOLVER_BACKEND_VEXCL)
294 BOOST_PP_SEQ_FOR_EACH(CALL_BLOCK_SOLVER, ~, ARCCORE_ALINA_BLOCK_SIZES)
295#endif
296
297 default:
298 precondition(false, "Unsupported block size");
299 break;
300 }
301
302 std::cout << prof << std::endl;
303}
Algebraic multigrid method.
Definition AMG.h:71
Convenience class that bundles together a preconditioner and an iterative solver.
Allows to use an AMG smoother as standalone preconditioner.
Definition Relaxation.h:144
Result of a solving.
Definition AlinaUtils.h:52
CPR preconditioner with Dynamic Row Sum modification.
Matrix class, to be used by user.
Vue d'un tableau d'éléments de type T.
Definition Span.h:801
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Alina::detail::empty_params params
Runtime-configurable wrappers around iterative solvers.