The ILU-preconditioned solver example..
This example depends on simple-solver.
Introduction
About the example
This example shows how to use incomplete factors generated via the ParILU algorithm to generate an incomplete factorization (ILU) preconditioner, how to specify the sparse triangular solves in the ILU preconditioner application, and how to generate an ILU-preconditioned solver and apply it to a specific problem.
The commented program
}
const auto executor_string = argc >= 2 ? argv[1] : "reference";
Figure out where to run the code
std::map<std::string, std::function<std::shared_ptr<gko::Executor>()>>
exec_map{
{"cuda",
[] {
}},
{"hip",
[] {
}},
{"dpcpp",
[] {
}},
{"reference", [] { return gko::ReferenceExecutor::create(); }}};
executor where Ginkgo will perform the computation
const auto exec = exec_map.at(executor_string)();
Read data
auto A =
gko::share(gko::read<mtx>(std::ifstream(
"data/A.mtx"), exec));
auto b = gko::read<vec>(std::ifstream("data/b.mtx"), exec);
auto x = gko::read<vec>(std::ifstream("data/x0.mtx"), exec);
Generate incomplete factors using ParILU
Generate concrete factorization for input matrix
auto par_ilu =
gko::share(par_ilu_fact->generate(A));
Generate an ILU preconditioner factory by setting lower and upper triangular solver - in this case the exact triangular solves
auto ilu_pre_factory =
false>::build()
.on(exec);
Use incomplete factors to generate ILU preconditioner
auto ilu_preconditioner =
gko::share(ilu_pre_factory->generate(par_ilu));
Use preconditioner inside GMRES solver factory Generating a solver factory tied to a specific preconditioner makes sense if there are several very similar systems to solve, and the same solver+preconditioner combination is expected to be effective.
const RealValueType reduction_factor{1e-7};
auto ilu_gmres_factory =
gmres::build()
.with_criteria(gko::stop::Iteration::build().with_max_iters(1000u),
.with_reduction_factor(reduction_factor))
.with_generated_preconditioner(ilu_preconditioner)
.on(exec);
Generate preconditioned solver for a specific target system
auto ilu_gmres = ilu_gmres_factory->generate(A);
Solve system
Print solution
std::cout << "Solution (x):\n";
Calculate residual
auto one = gko::initialize<vec>({1.0}, exec);
auto neg_one = gko::initialize<vec>({-1.0}, exec);
auto res = gko::initialize<real_vec>({0.0}, exec);
A->apply(one, x, neg_one, b);
b->compute_norm2(res);
std::cout << "Residual norm sqrt(r^T r):\n";
}
Results
This is the expected output:
Solution (x):
19 1
0.252218
0.108645
0.0662811
0.0630433
0.0384088
0.0396536
0.0402648
0.0338935
0.0193098
0.0234653
0.0211499
0.0196413
0.0199151
0.0181674
0.0162722
0.0150714
0.0107016
0.0121141
0.0123025
Residual norm sqrt(r^T r):
1 1
1.46249e-08
Comments about programming and debugging
The plain program
#include <ginkgo/ginkgo.hpp>
#include <cstdlib>
#include <fstream>
#include <iostream>
#include <map>
#include <string>
int main(int argc, char* argv[])
{
using ValueType = double;
using IndexType = int;
if (argc == 2 && (std::string(argv[1]) == "--help")) {
std::cerr << "Usage: " << argv[0] << " [executor]" << std::endl;
std::exit(-1);
}
const auto executor_string = argc >= 2 ? argv[1] : "reference";
std::map<std::string, std::function<std::shared_ptr<gko::Executor>()>>
exec_map{
{"cuda",
[] {
}},
{"hip",
[] {
}},
{"dpcpp",
[] {
}},
{"reference", [] { return gko::ReferenceExecutor::create(); }}};
const auto exec = exec_map.at(executor_string)();
auto A =
gko::share(gko::read<mtx>(std::ifstream(
"data/A.mtx"), exec));
auto b = gko::read<vec>(std::ifstream("data/b.mtx"), exec);
auto x = gko::read<vec>(std::ifstream("data/x0.mtx"), exec);
auto par_ilu_fact =
auto par_ilu =
gko::share(par_ilu_fact->generate(A));
auto ilu_pre_factory =
false>::build()
.on(exec);
auto ilu_preconditioner =
gko::share(ilu_pre_factory->generate(par_ilu));
const RealValueType reduction_factor{1e-7};
auto ilu_gmres_factory =
gmres::build()
.with_criteria(gko::stop::Iteration::build().with_max_iters(1000u),
.with_reduction_factor(reduction_factor))
.with_generated_preconditioner(ilu_preconditioner)
.on(exec);
auto ilu_gmres = ilu_gmres_factory->generate(A);
ilu_gmres->apply(b, x);
std::cout << "Solution (x):\n";
auto one = gko::initialize<vec>({1.0}, exec);
auto neg_one = gko::initialize<vec>({-1.0}, exec);
auto res = gko::initialize<real_vec>({0.0}, exec);
A->apply(one, x, neg_one, b);
b->compute_norm2(res);
std::cout << "Residual norm sqrt(r^T r):\n";
}