Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
NonScalarSearch.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
24#include "arccore/alina/BuiltinBackend.h"
25#include "arccore/alina/RelaxationRuntime.h"
26#include "arccore/alina/CoarseningRuntime.h"
27#include "arccore/alina/SolverRuntime.h"
28#include "arccore/alina/PreconditionerRuntime.h"
29#include "arccore/alina/DeflatedSolver.h"
30#include "arccore/alina/AMG.h"
31#include "arccore/alina/Adapters.h"
32#include "arccore/alina/IO.h"
33
34#include "arccore/alina/Profiler.h"
35
36using namespace Arcane;
37
38using Alina::precondition;
39
40//---------------------------------------------------------------------------
41int main(int argc, char* argv[])
42{
43 auto& prof = Alina::Profiler::globalProfiler();
44
45 namespace po = boost::program_options;
46 namespace io = Alina::IO;
47
48 using std::string;
49 using std::vector;
50
51 po::options_description desc("Options");
52
53 desc.add_options()("help,h", "Show this help.")("prm-file,P",
54 po::value<string>(),
55 "Parameter file in json format. ")(
56 "prm,p",
57 po::value<vector<string>>()->multitoken(),
58 "Parameters specified as name=value pairs. "
59 "May be provided multiple times. Examples:\n"
60 " -p solver.tol=1e-3\n"
61 " -p precond.coarse_enough=300")("matrix,A",
62 po::value<string>()->required(),
63 "System matrix in the MatrixMarket format.")(
64 "rhs,f",
65 po::value<string>(),
66 "The RHS vector in the MatrixMarket format. "
67 "When omitted, a vector of ones is used by default. "
68 "Should only be provided together with a system matrix. ")(
69 "scale,s",
70 po::bool_switch()->default_value(false),
71 "Scale the matrix so that the diagonal is unit. ")(
72 "null,N",
73 po::value<string>(),
74 "Starting null-vectors in the MatrixMarket format. ")(
75 "numvec,n",
76 po::value<int>()->default_value(3),
77 "The number of near nullspace vectors to search for. ")(
78 "binary,B",
79 po::bool_switch()->default_value(false),
80 "When specified, treat input files as binary instead of as MatrixMarket. "
81 "It is assumed the files were converted to binary format with mm2bin utility. ")(
82 "output,o",
83 po::value<string>(),
84 "Output the computed nullspace to the MatrixMarket file.");
85
86 po::positional_options_description p;
87 p.add("prm", -1);
88
89 po::variables_map vm;
90 po::store(po::command_line_parser(argc, argv).options(desc).positional(p).run(), vm);
91 po::notify(vm);
92
93 if (vm.count("help")) {
94 std::cout << desc << std::endl;
95 return 0;
96 }
97
98 for (int i = 0; i < argc; ++i) {
99 if (i)
100 std::cout << " ";
101 std::cout << argv[i];
102 }
103 std::cout << std::endl;
104
106 if (vm.count("prm-file")) {
107 prm.read_json(vm["prm-file"].as<string>());
108 }
109
110 if (vm.count("prm")) {
111 for (const string& v : vm["prm"].as<vector<string>>()) {
112 prm.putKeyValue(v);
113 }
114 }
115
116 ptrdiff_t rows, nv = 0, numvec = vm["numvec"].as<int>();
117 vector<ptrdiff_t> ptr, col;
118 vector<double> val, rhs;
119 std::list<std::vector<double>> Z;
120
121 {
122 auto t = prof.scoped_tic("read");
123
124 string Afile = vm["matrix"].as<string>();
125 bool binary = vm["binary"].as<bool>();
126
127 if (binary) {
128 io::read_crs(Afile, rows, ptr, col, val);
129 }
130 else {
131 ptrdiff_t cols;
132 std::tie(rows, cols) = io::mm_reader(Afile)(ptr, col, val);
133 precondition(rows == cols, "Non-square system matrix");
134 }
135
136 if (vm.count("rhs")) {
137 string bfile = vm["rhs"].as<string>();
138
139 ptrdiff_t n, m;
140
141 if (binary) {
142 io::read_dense(bfile, n, m, rhs);
143 }
144 else {
145 std::tie(n, m) = io::mm_reader(bfile)(rhs);
146 }
147
148 precondition(n == rows && m == 1, "The RHS vector has wrong size");
149 }
150 else {
151 rhs.resize(rows, 1.0);
152 }
153
154 if (vm.count("null")) {
155 string nfile = vm["null"].as<string>();
156
157 std::vector<double> null;
158 ptrdiff_t m;
159
160 if (binary) {
161 io::read_dense(nfile, m, nv, null);
162 }
163 else {
164 std::tie(m, nv) = io::mm_reader(nfile)(null);
165 }
166
167 precondition(m == rows, "Near null-space vectors have wrong size");
168
169 for (ptrdiff_t i = 0; i < nv; ++i) {
170 Z.emplace_back(rows);
171 for (ptrdiff_t j = 0; j < rows; ++j) {
172 Z.back()[j] = null[j * nv + i];
173 }
174 }
175 }
176 }
177
178 if (vm["scale"].as<bool>()) {
179 auto t = prof.scoped_tic("scaling");
180 std::vector<double> dia(rows, 1.0);
181
182 for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(rows); ++i) {
183 double d = 1.0;
184 for (ptrdiff_t j = ptr[i], e = ptr[i + 1]; j < e; ++j) {
185 if (col[j] == i) {
186 d = 1 / sqrt(val[j]);
187 }
188 }
189 if (!std::isnan(d))
190 dia[i] = d;
191 }
192
193 for (ptrdiff_t i = 0; i < static_cast<ptrdiff_t>(rows); ++i) {
194 rhs[i] *= dia[i];
195 for (ptrdiff_t j = ptr[i], e = ptr[i + 1]; j < e; ++j) {
196 val[j] *= dia[i] * dia[col[j]];
197 }
198 }
199 }
200
201 using Backend = Alina::BuiltinBackend<double>;
204
205 std::mt19937 rng;
206 std::uniform_real_distribution<double> rnd(-1, 1);
207 std::vector<double> x(rows), zero(rows, 0.0);
208
209 auto A = std::tie(rows, ptr, col, val);
210
211 prm.put("solver.ns_search", true);
212
213 prof.tic("search");
214 for (int k = nv; k < numvec; ++k) {
215 auto t = prof.scoped_tic(std::string("vector ") + std::to_string(k));
216 std::vector<double> N;
217
218 if (k) {
219 N.resize(k * rows);
220 int j = 0;
221 for (const auto& z : Z) {
222 for (ptrdiff_t i = 0; i < rows; ++i) {
223 N[i * k + j] = z[i];
224 }
225 ++j;
226 }
227
228 prm.put("precond.coarsening.nullspace.rows", rows);
229 prm.put("precond.coarsening.nullspace.cols", k);
230 prm.put("precond.coarsening.nullspace.B", N.data());
231 }
232
233 prof.tic("setup");
234 Solver S(A, prm);
235 prof.toc("setup");
236
237 std::cout << std::endl
238 << "-------------------------" << std::endl
239 << "-- Searching for vector " << k << std::endl
240 << "-------------------------" << std::endl
241 << S << std::endl;
242
243 for (auto& v : x)
244 v = rnd(rng);
245
246 prof.tic("solve");
247 Alina::SolverResult r = S(zero, x);
248 prof.toc("solve");
249
250 std::cout << "Iterations: " << r.nbIteration() << std::endl
251 << "Error: " << r.residual() << std::endl;
252
253 // Orthonormalize the new vector
254 for (const auto& z : Z) {
255 double c = Alina::backend::inner_product(x, z) / Alina::backend::inner_product(z, z);
256 Alina::backend::axpby(-c, z, 1, x);
257 }
258
259 double nx = sqrt(Alina::backend::inner_product(x, x));
260 for (auto& v : x)
261 v /= nx;
262 Z.push_back(x);
263 }
264 prof.toc("search");
265
266 // Solve the system using the near nullspace vectors:
267 std::vector<double> N(numvec * rows);
268 {
269 auto t = prof.scoped_tic("apply");
270
271 int j = 0;
272 for (const auto& z : Z) {
273 for (ptrdiff_t i = 0; i < rows; ++i) {
274 N[i * numvec + j] = z[i];
275 }
276 ++j;
277 }
278
279 prm.put("precond.coarsening.nullspace.rows", rows);
280 prm.put("precond.coarsening.nullspace.cols", numvec);
281 prm.put("precond.coarsening.nullspace.B", N.data());
282
283 prof.tic("setup");
284 Solver S(A, prm);
285 prof.toc("setup");
286
287 std::cout << std::endl
288 << "-------------------------" << std::endl
289 << "-- Solving the system " << std::endl
290 << "-------------------------" << std::endl
291 << S << std::endl;
292
293 Alina::backend::clear(x);
294
295 prof.tic("solve");
296 Alina::SolverResult r = S(rhs, x);
297 prof.toc("solve");
298
299 std::cout << "Iterations: " << r.nbIteration() << std::endl
300 << "Error: " << r.residual() << std::endl;
301 }
302
303 if (vm.count("output")) {
304 auto t = prof.scoped_tic("write");
305 Alina::IO::mm_write(vm["output"].as<string>(), N.data(), rows, numvec);
306 }
307
308 std::cout << prof << std::endl;
309}
Convenience class that bundles together a preconditioner and an iterative solver.
Result of a solving.
Definition AlinaUtils.h:52
__host__ __device__ double sqrt(double v)
Racine carrée de v.
Definition Math.h:135
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Runtime-configurable wrappers around iterative solvers.