Ginkgo  Generated from pipelines/1068515030 branch based on master. Ginkgo version 1.7.0
A numerical linear algebra library targeting many-core architectures
The ir-ilu-preconditioned-solver program

The IR-ILU preconditioned solver example..

This example depends on ilu-preconditioned-solver, iterative-refinement.

Table of contents
  1. Introduction
  2. The commented program
  1. Results
  2. The plain program

Introduction

About the example

This example shows how to combine iterative refinement with the adaptive precision block-Jacobi preconditioner in order to approximately solve the triangular systems occurring in ILU preconditioning. Using an adaptive precision block-Jacobi preconditioner matrix as inner solver for the iterative refinement method is equivalent to doing adaptive precision block-Jacobi relaxation in the triangular solves. This example roughly approximates the triangular solves with five adaptive precision block-Jacobi sweeps with a maximum block size of 16.

This example is motivated by "Multiprecision block-Jacobi for Iterative Triangular Solves" (Göbel, Anzt, Cojean, Flegar, Quintana-Ortí, Euro-Par 2020). The theory and a detailed analysis can be found there.

The commented program

#include <ginkgo/ginkgo.hpp>
#include <cstdlib>
#include <fstream>
#include <iostream>
#include <map>
#include <string>
int main(int argc, char* argv[])
{

Some shortcuts

using ValueType = double;
using RealValueType = gko::remove_complex<ValueType>;
using IndexType = int;

Print version information

std::cout << gko::version_info::get() << std::endl;

Figure out where to run the code

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";
const unsigned int sweeps = argc == 3 ? std::atoi(argv[2]) : 5u;
std::map<std::string, std::function<std::shared_ptr<gko::Executor>()>>
exec_map{
{"omp", [] { return gko::OmpExecutor::create(); }},
{"cuda",
[] {
}},
{"hip",
[] {
}},
{"dpcpp",
[] {
}},
{"reference", [] { return gko::ReferenceExecutor::create(); }}};

executor where Ginkgo will perform the computation

const auto exec = exec_map.at(executor_string)(); // throws if not valid

Read data

auto A = gko::share(gko::read<mtx>(std::ifstream("data/A.mtx"), exec));

Create RHS and initial guess as 1

gko::size_type num_rows = A->get_size()[0];
auto host_x = vec::create(exec->get_master(), gko::dim<2>(num_rows, 1));
for (gko::size_type i = 0; i < num_rows; i++) {
host_x->at(i, 0) = 1.;
}
auto x = gko::clone(exec, host_x);
auto b = gko::clone(exec, host_x);
auto clone_x = gko::clone(exec, x);

Generate incomplete factors using ParILU

Generate concrete factorization for input matrix

auto par_ilu = gko::share(par_ilu_fact->generate(A));

Generate an iterative refinement factory to be used as a triangular solver in the preconditioner application. The generated method is equivalent to doing five block-Jacobi sweeps with a maximum block size of 16.

auto bj_factory = gko::share(
bj::build()
.with_max_block_size(16u)
.with_storage_optimization(gko::precision_reduction::autodetect())
.on(exec));
auto trisolve_factory =
ir::build()
.with_solver(bj_factory)
.with_criteria(gko::stop::Iteration::build().with_max_iters(sweeps))
.on(exec);

Generate an ILU preconditioner factory by setting lower and upper triangular solver - in this case the previously defined iterative refinement method.

.with_l_solver(gko::clone(trisolve_factory))
.with_u_solver(gko::clone(trisolve_factory))
.on(exec);

Use incomplete factors to generate ILU preconditioner

auto ilu_preconditioner = gko::share(ilu_pre_factory->generate(par_ilu));

Create stopping criteria for Gmres

const RealValueType reduction_factor{1e-12};
auto iter_stop = gko::share(
gko::stop::Iteration::build().with_max_iters(1000u).on(exec));
.with_reduction_factor(reduction_factor)
.on(exec));

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.

auto ilu_gmres_factory =
gmres::build()
.with_criteria(iter_stop, tol_stop)
.with_generated_preconditioner(ilu_preconditioner)
.on(exec);

Generate preconditioned solver for a specific target system

auto ilu_gmres = ilu_gmres_factory->generate(A);

Add logger

std::shared_ptr<const gko::log::Convergence<ValueType>> logger =
ilu_gmres->add_logger(logger);

Warmup run

ilu_gmres->apply(b, x);

Solve system 100 times and take the average time.

std::chrono::nanoseconds time(0);
for (int i = 0; i < 100; i++) {
x->copy_from(clone_x);
auto tic = std::chrono::high_resolution_clock::now();
ilu_gmres->apply(b, x);
auto toc = std::chrono::high_resolution_clock::now();
time += std::chrono::duration_cast<std::chrono::nanoseconds>(toc - tic);
}
std::cout << "Using " << sweeps << " block-Jacobi sweeps.\n";

Print solution

std::cout << "Solution (x):\n";
write(std::cout, x);

Get residual

auto res = gko::as<vec>(logger->get_residual_norm());
std::cout << "GMRES iteration count: " << logger->get_num_iterations()
<< "\n";
std::cout << "GMRES execution time [ms]: "
<< static_cast<double>(time.count()) / 100000000.0 << "\n";
std::cout << "Residual norm sqrt(r^T r):\n";
write(std::cout, res);
}

Results

This is the expected output:

Using 5 block-Jacobi sweeps.
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
GMRES iteration count: 8
GMRES execution time [ms]: 0.377673
Residual norm sqrt(r^T r):
%%MatrixMarket matrix array real general
1 1
1.65303e-12

Comments about programming and debugging

The plain program

/*******************************<GINKGO LICENSE>******************************
Copyright (c) 2017-2023, the Ginkgo authors
All rights reserved.
Redistribution and use in source and binary forms, with or without
modification, are permitted provided that the following conditions
are met:
1. Redistributions of source code must retain the above copyright
notice, this list of conditions and the following disclaimer.
2. Redistributions in binary form must reproduce the above copyright
notice, this list of conditions and the following disclaimer in the
documentation and/or other materials provided with the distribution.
3. Neither the name of the copyright holder nor the names of its
contributors may be used to endorse or promote products derived from
this software without specific prior written permission.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS
IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A
PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
(INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
******************************<GINKGO LICENSE>*******************************/
#include <ginkgo/ginkgo.hpp>
#include <cstdlib>
#include <fstream>
#include <iostream>
#include <map>
#include <string>
int main(int argc, char* argv[])
{
using ValueType = double;
using RealValueType = gko::remove_complex<ValueType>;
using IndexType = int;
std::cout << gko::version_info::get() << std::endl;
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";
const unsigned int sweeps = argc == 3 ? std::atoi(argv[2]) : 5u;
std::map<std::string, std::function<std::shared_ptr<gko::Executor>()>>
exec_map{
{"omp", [] { return gko::OmpExecutor::create(); }},
{"cuda",
[] {
}},
{"hip",
[] {
}},
{"dpcpp",
[] {
}},
{"reference", [] { return gko::ReferenceExecutor::create(); }}};
const auto exec = exec_map.at(executor_string)(); // throws if not valid
auto A = gko::share(gko::read<mtx>(std::ifstream("data/A.mtx"), exec));
gko::size_type num_rows = A->get_size()[0];
auto host_x = vec::create(exec->get_master(), gko::dim<2>(num_rows, 1));
for (gko::size_type i = 0; i < num_rows; i++) {
host_x->at(i, 0) = 1.;
}
auto x = gko::clone(exec, host_x);
auto b = gko::clone(exec, host_x);
auto clone_x = gko::clone(exec, x);
auto par_ilu_fact =
auto par_ilu = gko::share(par_ilu_fact->generate(A));
auto bj_factory = gko::share(
bj::build()
.with_max_block_size(16u)
.with_storage_optimization(gko::precision_reduction::autodetect())
.on(exec));
auto trisolve_factory =
ir::build()
.with_solver(bj_factory)
.with_criteria(gko::stop::Iteration::build().with_max_iters(sweeps))
.on(exec);
auto ilu_pre_factory = gko::preconditioner::Ilu<ir, ir>::build()
.with_l_solver(gko::clone(trisolve_factory))
.with_u_solver(gko::clone(trisolve_factory))
.on(exec);
auto ilu_preconditioner = gko::share(ilu_pre_factory->generate(par_ilu));
const RealValueType reduction_factor{1e-12};
auto iter_stop = gko::share(
gko::stop::Iteration::build().with_max_iters(1000u).on(exec));
.with_reduction_factor(reduction_factor)
.on(exec));
auto ilu_gmres_factory =
gmres::build()
.with_criteria(iter_stop, tol_stop)
.with_generated_preconditioner(ilu_preconditioner)
.on(exec);
auto ilu_gmres = ilu_gmres_factory->generate(A);
std::shared_ptr<const gko::log::Convergence<ValueType>> logger =
ilu_gmres->add_logger(logger);
ilu_gmres->apply(b, x);
std::chrono::nanoseconds time(0);
for (int i = 0; i < 100; i++) {
x->copy_from(clone_x);
auto tic = std::chrono::high_resolution_clock::now();
ilu_gmres->apply(b, x);
auto toc = std::chrono::high_resolution_clock::now();
time += std::chrono::duration_cast<std::chrono::nanoseconds>(toc - tic);
}
std::cout << "Using " << sweeps << " block-Jacobi sweeps.\n";
std::cout << "Solution (x):\n";
write(std::cout, x);
auto res = gko::as<vec>(logger->get_residual_norm());
std::cout << "GMRES iteration count: " << logger->get_num_iterations()
<< "\n";
std::cout << "GMRES execution time [ms]: "
<< static_cast<double>(time.count()) / 100000000.0 << "\n";
std::cout << "Residual norm sqrt(r^T r):\n";
write(std::cout, res);
}
gko::matrix::Csr
CSR is a matrix format which stores only the nonzero coefficients by compressing each row of the matr...
Definition: matrix.hpp:54
gko::log::Convergence::create
static std::unique_ptr< Convergence > create(std::shared_ptr< const Executor >, const mask_type &enabled_events=Logger::criterion_events_mask|Logger::iteration_complete_mask)
Creates a convergence logger.
Definition: convergence.hpp:106
gko::layout_type::array
The matrix should be written as dense matrix in column-major order.
gko::matrix::Dense
Dense is a matrix format which explicitly stores all values of the matrix.
Definition: dense_cache.hpp:48
gko::precision_reduction::autodetect
constexpr static precision_reduction autodetect() noexcept
Returns a special encoding which instructs the algorithm to automatically detect the best precision.
Definition: types.hpp:347
gko::solver::Ir
Iterative refinement (IR) is an iterative method that uses another coarse method to approximate the e...
Definition: ir.hpp:108
gko::size_type
std::size_t size_type
Integral type used for allocation quantities.
Definition: types.hpp:120
gko::solver::Gmres
GMRES or the generalized minimal residual method is an iterative type Krylov subspace method which is...
Definition: gmres.hpp:76
gko::clone
detail::cloned_type< Pointer > clone(const Pointer &p)
Creates a unique clone of the object pointed to by p.
Definition: utils_helper.hpp:203
gko::HipExecutor::create
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.
gko::version_info::get
static const version_info & get()
Returns an instance of version_info.
Definition: version.hpp:168
gko::preconditioner::Jacobi
A block-Jacobi preconditioner is a block-diagonal linear operator, obtained by inverting the diagonal...
Definition: jacobi.hpp:213
gko::stop::ResidualNorm
The ResidualNorm class is a stopping criterion which stops the iteration process when the actual resi...
Definition: residual_norm.hpp:138
gko::dim< 2 >
gko::write
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:324
gko::share
detail::shared_type< OwningPointer > share(OwningPointer &&p)
Marks the object pointed to by p as shared.
Definition: utils_helper.hpp:254
gko::factorization::ParIlu
ParILU is an incomplete LU factorization which is computed in parallel.
Definition: par_ilu.hpp:97
gko::CudaExecutor::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.
gko::OmpExecutor::create
static std::shared_ptr< OmpExecutor > create(std::shared_ptr< CpuAllocatorBase > alloc=std::make_shared< CpuAllocator >())
Creates a new OmpExecutor.
Definition: executor.hpp:1373
gko::remove_complex
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:354
gko::DpcppExecutor::create
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.
gko::real
constexpr auto real(const T &x)
Returns the real part of the object.
Definition: math.hpp:1041
gko::preconditioner::Ilu
The Incomplete LU (ILU) preconditioner solves the equation for a given lower triangular matrix L,...
Definition: ilu.hpp:113