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 <<
"Printing solver infos\n";
105 std::cout << solve << std::endl;
107 solve.prm.get(ptree);
108 std::cout <<
"SOLVER PARAMS: " << ptree <<
"\n";
110 auto f_b = Backend::copy_vector(rhs, bprm);
111 auto x_b = Backend::copy_vector(x, bprm);
114 info = solve(*f_b, *x_b);
127 auto& prof = Alina::Profiler::globalProfiler();
128 namespace po = boost::program_options;
129 namespace io = Alina::IO;
134 po::options_description desc(
"Options");
136 desc.add_options()(
"help,h",
"Show this help.");
137 desc.add_options()(
"prm-file,P", po::value<string>(),
138 "Parameter file in json format. ");
139 desc.add_options()(
"prm,p", po::value<vector<string>>()->multitoken(),
140 "Parameters specified as name=value pairs. "
141 "May be provided multiple times. Examples:\n"
142 " -p solver.tol=1e-3\n"
143 " -p precond.coarse_enough=300");
145 desc.add_options()(
"size,n",
146 po::value<int>()->default_value(32),
147 "The size of the Poisson problem to solve when no system matrix is given. "
148 "Specified as number of grid nodes along each dimension of a unit cube. "
149 "The resulting system will have n*n*n unknowns. ");
151 desc.add_options()(
"anisotropy,a",
152 po::value<double>()->default_value(1.0),
153 "The anisotropy value for the generated Poisson value. "
154 "Used to determine problem scaling along X, Y, and Z axes: "
155 "hy = hx * a, hz = hy * a.");
157 po::positional_options_description p;
160 po::variables_map vm;
161 po::store(po::command_line_parser(argc, argv).options(desc).positional(p).run(), vm);
164 if (vm.count(
"help")) {
169 for (
int i = 0; i < argc; ++i) {
172 std::cout << argv[i];
174 std::cout << std::endl;
177 if (vm.count(
"prm-file")) {
178 prm.read_json(vm[
"prm-file"].as<string>());
181 if (vm.count(
"prm")) {
182 for (
const string& v : vm[
"prm"].as<vector<string>>()) {
187 size_t rows = vm[
"size"].as<
int>();
188 vector<ptrdiff_t> ptr, col;
189 vector<double> val, rhs, null, x;
190 std::cout <<
"ROWS=" << rows <<
"\n";
192 auto t = prof.scoped_tic(
"assembling");
193 rows = sample_problem(rows, val, col, ptr, rhs, vm[
"anisotropy"].as<double>());
199 bool do_hypre =
false;
201 do_hypre = v.value();
205 _doHypreSolver(rows, ptr, col, val, rhs, x, argc, argv);
209 solver_result = solve(prm, rows, ptr, col, val, rhs, x);
211 if (vm.count(
"output")) {
212 auto t = prof.scoped_tic(
"write");
213 Alina::IO::mm_write(vm[
"output"].as<string>(), &x[0], x.size());
216 std::cout <<
"Iterations: " << solver_result.nbIteration() << std::endl
217 <<
"Error: " << solver_result.residual() << std::endl
218 << prof << std::endl;
225int main(
int argc,
char* argv[])
227 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