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(); }}};
static std::shared_ptr< CudaExecutor > create(int device_id, std::shared_ptr< Executor > master, bool device_reset, allocation_mode alloc_mode=default_cuda_alloc_mode, CUstream_st *stream=nullptr)
Creates a new CudaExecutor.
static std::shared_ptr< DpcppExecutor > create(int device_id, std::shared_ptr< Executor > master, std::string device_type="all", dpcpp_queue_property property=dpcpp_queue_property::in_order)
Creates a new DpcppExecutor.
static std::shared_ptr< HipExecutor > create(int device_id, std::shared_ptr< Executor > master, bool device_reset, allocation_mode alloc_mode=default_hip_alloc_mode, CUstream_st *stream=nullptr)
Creates a new HipExecutor.
static std::shared_ptr< OmpExecutor > create(std::shared_ptr< CpuAllocatorBase > alloc=std::make_shared< CpuAllocator >())
Creates a new OmpExecutor.
Definition executor.hpp:1396
executor where Ginkgo will perform the computation
const auto exec = exec_map.at(executor_string)();
Read data
std::unique_ptr< MatrixType > read(StreamType &&is, MatrixArgs &&... args)
Reads a matrix stored in matrix market format from an input stream.
Definition mtx_io.hpp:159
detail::shared_type< OwningPointer > share(OwningPointer &&p)
Marks the object pointed to by p as shared.
Definition utils_helper.hpp:224
Generate incomplete factors using ParILU
auto par_ilu_fact =
ParILU is an incomplete LU factorization which is computed in parallel.
Definition par_ilu.hpp:70
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);
The Incomplete LU (ILU) preconditioner solves the equation for a given lower triangular matrix L,...
Definition ilu.hpp:124
UpperTrs is the triangular solver which solves the system U x = b, when U is an upper triangular matr...
Definition triangular.hpp:236
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);
The ResidualNorm class is a stopping criterion which stops the iteration process when the actual resi...
Definition residual_norm.hpp:113
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";
write(std::cout, x);
Calculate residual
A->apply(one, x, neg_one, b);
b->compute_norm2(res);
std::cout << "Residual norm sqrt(r^T r):\n";
}
std::unique_ptr< Matrix > initialize(size_type stride, std::initializer_list< typename Matrix::value_type > vals, std::shared_ptr< const Executor > exec, TArgs &&... create_args)
Creates and initializes a column-vector.
Definition dense.hpp:1565
void write(StreamType &&os, MatrixPtrType &&matrix, layout_type layout=detail::mtx_io_traits< std::remove_cv_t< detail::pointee< MatrixPtrType > > >::default_layout)
Writes a matrix into an output stream in matrix market format.
Definition mtx_io.hpp:295
Results
This is the expected output:
Solution (x):
%%MatrixMarket matrix array real general
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):
%%MatrixMarket matrix array real general
1 1
1.46249e-08
Comments about programming and debugging
The plain program
#include <cstdlib>
#include <fstream>
#include <iostream>
#include <map>
#include <string>
#include <ginkgo/ginkgo.hpp>
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 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";
A->apply(one, x, neg_one, b);
b->compute_norm2(res);
std::cout << "Residual norm sqrt(r^T r):\n";
}
CSR is a matrix format which stores only the nonzero coefficients by compressing each row of the matr...
Definition sparsity_csr.hpp:21
Dense is a matrix format which explicitly stores all values of the matrix.
Definition sparsity_csr.hpp:25
GMRES or the generalized minimal residual method is an iterative type Krylov subspace method which is...
Definition gmres.hpp:75
static const version_info & get()
Returns an instance of version_info.
Definition version.hpp:139
constexpr T one()
Returns the multiplicative identity for T.
Definition math.hpp:630
typename detail::remove_complex_s< T >::type remove_complex
Obtain the type which removed the complex of complex/scalar type or the template parameter of class b...
Definition math.hpp:260