22#include <boost/program_options.hpp>
24#include "arccore/alina/BuiltinBackend.h"
25#include "arccore/alina/StaticMatrix.h"
26#include "arccore/alina/Adapters.h"
28#include "arcane/utils/PlatformUtils.h"
29#include "arcane/utils/String.h"
30#include "arcane/utils/Convert.h"
32#include "arccore/alina/Relaxation.h"
33#include "arccore/alina/Coarsening.h"
34#include "arccore/alina/BiCGStabSolver.h"
35#include "arccore/alina/PreconditionedSolver.h"
36#include "arccore/alina/AMG.h"
37#include "arccore/alina/Adapters.h"
38#include "arccore/alina/IO.h"
40#include "arccore/alina/SolverRuntime.h"
41#include "arccore/alina/PreconditionerRuntime.h"
43#include "arccore/alina/Profiler.h"
45#include "AlinaSamplesCommon.h"
46#include "arccore/trace/ITraceMng.h"
48#include "SampleProblemCommon.h"
51using Alina::precondition;
63_doHypreSolver(
int nb_row,
64 std::vector<ptrdiff_t>
const& ptr,
65 std::vector<ptrdiff_t>
const& col,
66 std::vector<double>
const& val,
67 std::vector<double>
const& rhs,
68 std::vector<double>& x,
69 int argc,
char* argv[]);
77 std::vector<ptrdiff_t>
const& ptr,
78 std::vector<ptrdiff_t>
const& col,
79 std::vector<double>
const& val,
80 std::vector<double>
const& rhs,
81 std::vector<double>& x)
83 std::cout <<
"Using scalar solve ptr_size=" <<
sizeof(ptrdiff_t)
84 <<
" ptr_type_size=" <<
sizeof(Backend::ptr_type)
85 <<
" col_type_size=" <<
sizeof(Backend::col_type)
86 <<
" value_type_size=" <<
sizeof(Backend::value_type)
88 auto& prof = Alina::Profiler::globalProfiler();
101 Solver solve(std::tie(rows, ptr, col, val), prm, bprm);
104 std::cout << solve << std::endl;
106 auto f_b = Backend::copy_vector(rhs, bprm);
107 auto x_b = Backend::copy_vector(x, bprm);
110 info = solve(*f_b, *x_b);
123 auto& prof = Alina::Profiler::globalProfiler();
124 namespace po = boost::program_options;
125 namespace io = Alina::IO;
130 po::options_description desc(
"Options");
132 desc.add_options()(
"help,h",
"Show this help.");
133 desc.add_options()(
"prm-file,P", po::value<string>(),
134 "Parameter file in json format. ");
135 desc.add_options()(
"prm,p", po::value<vector<string>>()->multitoken(),
136 "Parameters specified as name=value pairs. "
137 "May be provided multiple times. Examples:\n"
138 " -p solver.tol=1e-3\n"
139 " -p precond.coarse_enough=300");
141 desc.add_options()(
"size,n",
142 po::value<int>()->default_value(32),
143 "The size of the Poisson problem to solve when no system matrix is given. "
144 "Specified as number of grid nodes along each dimension of a unit cube. "
145 "The resulting system will have n*n*n unknowns. ");
147 desc.add_options()(
"anisotropy,a",
148 po::value<double>()->default_value(1.0),
149 "The anisotropy value for the generated Poisson value. "
150 "Used to determine problem scaling along X, Y, and Z axes: "
151 "hy = hx * a, hz = hy * a.");
153 po::positional_options_description p;
156 po::variables_map vm;
157 po::store(po::command_line_parser(argc, argv).options(desc).positional(p).run(), vm);
160 if (vm.count(
"help")) {
165 for (
int i = 0; i < argc; ++i) {
168 std::cout << argv[i];
170 std::cout << std::endl;
173 if (vm.count(
"prm-file")) {
174 prm.read_json(vm[
"prm-file"].as<string>());
177 if (vm.count(
"prm")) {
178 for (
const string& v : vm[
"prm"].as<vector<string>>()) {
183 size_t rows = vm[
"size"].as<
int>();
184 vector<ptrdiff_t> ptr, col;
185 vector<double> val, rhs, null, x;
186 std::cout <<
"ROWS=" << rows <<
"\n";
188 auto t = prof.scoped_tic(
"assembling");
189 rows = sample_problem(rows, val, col, ptr, rhs, vm[
"anisotropy"].as<double>());
195 bool do_hypre =
false;
197 do_hypre = v.value();
201 _doHypreSolver(rows, ptr, col, val, rhs, x, argc, argv);
204 solver_result = solve(prm, rows, ptr, col, val, rhs, x);
206 if (vm.count(
"output")) {
207 auto t = prof.scoped_tic(
"write");
208 Alina::IO::mm_write(vm[
"output"].as<string>(), &x[0], x.size());
211 std::cout <<
"Iterations: " << solver_result.nbIteration() << std::endl
212 <<
"Error: " << solver_result.residual() << std::endl
213 << prof << std::endl;
220int main(
int argc,
char* argv[])
222 return Arcane::Alina::SampleMainContext::execMain(main2, argc, argv);
Algebraic multigrid method.
Convenience class that bundles together a preconditioner and an iterative solver.
Classe template pour convertir un type.
Interface du gestionnaire de traces.
virtual TraceMessage info()=0
Flot pour un message d'information.
Chaîne de caractères unicode.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
Alina::detail::empty_params params