Ginkgo  Generated from pipelines/1330831941 branch based on master. Ginkgo version 1.8.0
A numerical linear algebra library targeting many-core architectures
The preconditioner-export program

The preconditioner export example..

This example depends on simple-solver.

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

Introduction

About the example

This example shows how to explicitly generate and store preconditioners for a given matrix. It can also be used to inspect and debug the preconditioner generation.

The commented program

std::string name)
{
std::ofstream stream{name};
std::cerr << "Writing " << name << std::endl;
}
template <typename Function>
auto try_generate(Function fun) -> decltype(fun())
{
decltype(fun()) result;
try {
result = fun();
} catch (const gko::Error& err) {
std::cerr << "Error: " << err.what() << '\n';
std::exit(-1);
}
return result;
}
int main(int argc, char* argv[])
{

print usage message

if (argc < 2 || executors.find(argv[1]) == executors.end()) {
std::cerr << "Usage: executable"
<< " <reference|omp|cuda|hip|dpcpp> [<matrix-file>] "
"[<jacobi|ilu|parilu|parilut|ilu-isai|parilu-isai|parilut-"
"isai] [<preconditioner args>]\n";
std::cerr << "Jacobi parameters: [<max-block-size>] [<accuracy>] "
"[<storage-optimization:auto|0|1|2>]\n";
std::cerr << "ParILU parameters: [<iteration-count>]\n";
std::cerr
<< "ParILUT parameters: [<iteration-count>] [<fill-in-limit>]\n";
std::cerr << "ILU-ISAI parameters: [<sparsity-power>]\n";
std::cerr << "ParILU-ISAI parameters: [<iteration-count>] "
"[<sparsity-power>]\n";
std::cerr << "ParILUT-ISAI parameters: [<iteration-count>] "
"[<fill-in-limit>] [<sparsity-power>]\n";
return -1;
}

generate executor based on first argument

auto exec = try_generate([&] { return executors.at(argv[1])(); });

set matrix and preconditioner name with default values

std::string matrix = argc < 3 ? "data/A.mtx" : argv[2];
std::string precond = argc < 4 ? "jacobi" : argv[3];

load matrix file into Csr format

auto mtx = gko::share(try_generate([&] {
std::ifstream mtx_stream{matrix};
if (!mtx_stream) {
throw GKO_STREAM_ERROR("Unable to open matrix file");
}
std::cerr << "Reading " << matrix << std::endl;
return gko::read<gko::matrix::Csr<>>(mtx_stream, exec);
}));

concatenate remaining arguments for filename

std::string output_suffix;
for (auto i = 4; i < argc; ++i) {
output_suffix = output_suffix + "-" + argv[i];
}

handle different preconditioners

if (precond == "jacobi") {

jacobi: max_block_size, accuracy, storage_optimization

auto factory_parameter = gko::preconditioner::Jacobi<>::build();
if (argc >= 5) {
factory_parameter.with_max_block_size(
static_cast<gko::uint32>(std::stoi(argv[4])));
}
if (argc >= 6) {
factory_parameter.with_accuracy(std::stod(argv[5]));
}
if (argc >= 7) {
factory_parameter.with_storage_optimization(
std::string{argv[6]} == "auto"
: gko::precision_reduction(0, std::stoi(argv[6])));
}
auto factory = factory_parameter.on(exec);
auto jacobi = try_generate([&] { return factory->generate(mtx); });
output(jacobi, matrix + ".jacobi" + output_suffix);
} else if (precond == "ilu") {

ilu: no parameters

auto ilu = gko::as<gko::Composition<>>(try_generate([&] {
return gko::factorization::Ilu<>::build().on(exec)->generate(mtx);
}));
output(gko::as<gko::matrix::Csr<>>(ilu->get_operators()[0]),
matrix + ".ilu-l");
output(gko::as<gko::matrix::Csr<>>(ilu->get_operators()[1]),
matrix + ".ilu-u");
} else if (precond == "parilu") {

parilu: iterations

auto factory_parameter = gko::factorization::ParIlu<>::build();
if (argc >= 5) {
factory_parameter.with_iterations(
static_cast<gko::size_type>(std::stoi(argv[4])));
}
auto factory = factory_parameter.on(exec);
auto ilu = gko::as<gko::Composition<>>(
try_generate([&] { return factory->generate(mtx); }));
output(gko::as<gko::matrix::Csr<>>(ilu->get_operators()[0]),
matrix + ".parilu" + output_suffix + "-l");
output(gko::as<gko::matrix::Csr<>>(ilu->get_operators()[1]),
matrix + ".parilu" + output_suffix + "-u");
} else if (precond == "parilut") {

parilut: iterations, fill-in limit

auto factory_parameter = gko::factorization::ParIlut<>::build();
if (argc >= 5) {
factory_parameter.with_iterations(
static_cast<gko::size_type>(std::stoi(argv[4])));
}
if (argc >= 6) {
factory_parameter.with_fill_in_limit(std::stod(argv[5]));
}
auto factory = factory_parameter.on(exec);
auto ilut = gko::as<gko::Composition<>>(
try_generate([&] { return factory->generate(mtx); }));
output(gko::as<gko::matrix::Csr<>>(ilut->get_operators()[0]),
matrix + ".parilut" + output_suffix + "-l");
output(gko::as<gko::matrix::Csr<>>(ilut->get_operators()[1]),
matrix + ".parilut" + output_suffix + "-u");
} else if (precond == "ilu-isai") {

ilu-isai: sparsity power

auto fact_factory =
int sparsity_power = 1;
if (argc >= 5) {
sparsity_power = std::stoi(argv[4]);
}
auto factory =
.with_factorization(fact_factory)
.with_sparsity_power(sparsity_power))
.with_sparsity_power(sparsity_power))
.on(exec);
auto ilu_isai = try_generate([&] { return factory->generate(mtx); });
output(ilu_isai->get_l_solver()->get_approximate_inverse(),
matrix + ".ilu-isai" + output_suffix + "-l");
output(ilu_isai->get_u_solver()->get_approximate_inverse(),
matrix + ".ilu-isai" + output_suffix + "-u");
} else if (precond == "parilu-isai") {

parilu-isai: iterations, sparsity power

auto fact_parameter = gko::factorization::ParIlu<>::build();
int sparsity_power = 1;
if (argc >= 5) {
fact_parameter.with_iterations(
static_cast<gko::size_type>(std::stoi(argv[4])));
}
if (argc >= 6) {
sparsity_power = std::stoi(argv[5]);
}
auto fact_factory = gko::share(fact_parameter.on(exec));
auto factory =
.with_factorization(fact_factory)
.with_sparsity_power(sparsity_power))
.with_sparsity_power(sparsity_power))
.on(exec);
auto ilu_isai = try_generate([&] { return factory->generate(mtx); });
output(ilu_isai->get_l_solver()->get_approximate_inverse(),
matrix + ".parilu-isai" + output_suffix + "-l");
output(ilu_isai->get_u_solver()->get_approximate_inverse(),
matrix + ".parilu-isai" + output_suffix + "-u");
} else if (precond == "parilut-isai") {

parilut-isai: iterations, fill-in limit, sparsity power

auto fact_parameter = gko::factorization::ParIlut<>::build();
int sparsity_power = 1;
if (argc >= 5) {
fact_parameter.with_iterations(
static_cast<gko::size_type>(std::stoi(argv[4])));
}
if (argc >= 6) {
fact_parameter.with_fill_in_limit(std::stod(argv[5]));
}
if (argc >= 7) {
sparsity_power = std::stoi(argv[6]);
}
auto fact_factory = gko::share(fact_parameter.on(exec));
auto factory =
.with_factorization(fact_factory)
.with_sparsity_power(sparsity_power))
.with_sparsity_power(sparsity_power))
.on(exec);
auto ilu_isai = try_generate([&] { return factory->generate(mtx); });
output(ilu_isai->get_l_solver()->get_approximate_inverse(),
matrix + ".parilut-isai" + output_suffix + "-l");
output(ilu_isai->get_u_solver()->get_approximate_inverse(),
matrix + ".parilut-isai" + output_suffix + "-u");
}
}

Results

This is the expected output:

Usage: ./preconditioner-export <reference|omp|cuda|hip|dpcpp> [<matrix-file>] [<jacobi|ilu|parilu|parilut|ilu-isai|parilu-isai|parilut-isai] [<preconditioner args>]
Jacobi parameters: [<max-block-size>] [<accuracy>] [<storage-optimization:auto|0|1|2>]
ParILU parameters: [<iteration-count>]
ParILUT parameters: [<iteration-count>] [<fill-in-limit>]
ILU-ISAI parameters: [<sparsity-power>]
ParILU-ISAI parameters: [<iteration-count>] [<sparsity-power>]
ParILUT-ISAI parameters: [<iteration-count>] [<fill-in-limit>] [<sparsity-power>]

When specifying an executor:

Reading data/A.mtx
Writing data/A.mtx.jacobi

Comments about programming and debugging

The plain program

#include <ginkgo/ginkgo.hpp>
#include <fstream>
#include <functional>
#include <iostream>
#include <map>
#include <memory>
#include <string>
const std::map<std::string, std::function<std::shared_ptr<gko::Executor>()>>
executors{{"reference", [] { return gko::ReferenceExecutor::create(); }},
{"omp", [] { return gko::OmpExecutor::create(); }},
{"cuda",
[] {
0, gko::ReferenceExecutor::create());
}},
{"hip",
[] {
0, gko::ReferenceExecutor::create());
}},
{"dpcpp", [] {
0, gko::ReferenceExecutor::create());
}}};
std::string name)
{
std::ofstream stream{name};
std::cerr << "Writing " << name << std::endl;
}
template <typename Function>
auto try_generate(Function fun) -> decltype(fun())
{
decltype(fun()) result;
try {
result = fun();
} catch (const gko::Error& err) {
std::cerr << "Error: " << err.what() << '\n';
std::exit(-1);
}
return result;
}
int main(int argc, char* argv[])
{
if (argc < 2 || executors.find(argv[1]) == executors.end()) {
std::cerr << "Usage: executable"
<< " <reference|omp|cuda|hip|dpcpp> [<matrix-file>] "
"[<jacobi|ilu|parilu|parilut|ilu-isai|parilu-isai|parilut-"
"isai] [<preconditioner args>]\n";
std::cerr << "Jacobi parameters: [<max-block-size>] [<accuracy>] "
"[<storage-optimization:auto|0|1|2>]\n";
std::cerr << "ParILU parameters: [<iteration-count>]\n";
std::cerr
<< "ParILUT parameters: [<iteration-count>] [<fill-in-limit>]\n";
std::cerr << "ILU-ISAI parameters: [<sparsity-power>]\n";
std::cerr << "ParILU-ISAI parameters: [<iteration-count>] "
"[<sparsity-power>]\n";
std::cerr << "ParILUT-ISAI parameters: [<iteration-count>] "
"[<fill-in-limit>] [<sparsity-power>]\n";
return -1;
}
auto exec = try_generate([&] { return executors.at(argv[1])(); });
std::string matrix = argc < 3 ? "data/A.mtx" : argv[2];
std::string precond = argc < 4 ? "jacobi" : argv[3];
auto mtx = gko::share(try_generate([&] {
std::ifstream mtx_stream{matrix};
if (!mtx_stream) {
throw GKO_STREAM_ERROR("Unable to open matrix file");
}
std::cerr << "Reading " << matrix << std::endl;
return gko::read<gko::matrix::Csr<>>(mtx_stream, exec);
}));
std::string output_suffix;
for (auto i = 4; i < argc; ++i) {
output_suffix = output_suffix + "-" + argv[i];
}
if (precond == "jacobi") {
auto factory_parameter = gko::preconditioner::Jacobi<>::build();
if (argc >= 5) {
factory_parameter.with_max_block_size(
static_cast<gko::uint32>(std::stoi(argv[4])));
}
if (argc >= 6) {
factory_parameter.with_accuracy(std::stod(argv[5]));
}
if (argc >= 7) {
factory_parameter.with_storage_optimization(
std::string{argv[6]} == "auto"
: gko::precision_reduction(0, std::stoi(argv[6])));
}
auto factory = factory_parameter.on(exec);
auto jacobi = try_generate([&] { return factory->generate(mtx); });
output(jacobi, matrix + ".jacobi" + output_suffix);
} else if (precond == "ilu") {
auto ilu = gko::as<gko::Composition<>>(try_generate([&] {
return gko::factorization::Ilu<>::build().on(exec)->generate(mtx);
}));
output(gko::as<gko::matrix::Csr<>>(ilu->get_operators()[0]),
matrix + ".ilu-l");
output(gko::as<gko::matrix::Csr<>>(ilu->get_operators()[1]),
matrix + ".ilu-u");
} else if (precond == "parilu") {
auto factory_parameter = gko::factorization::ParIlu<>::build();
if (argc >= 5) {
factory_parameter.with_iterations(
static_cast<gko::size_type>(std::stoi(argv[4])));
}
auto factory = factory_parameter.on(exec);
auto ilu = gko::as<gko::Composition<>>(
try_generate([&] { return factory->generate(mtx); }));
output(gko::as<gko::matrix::Csr<>>(ilu->get_operators()[0]),
matrix + ".parilu" + output_suffix + "-l");
output(gko::as<gko::matrix::Csr<>>(ilu->get_operators()[1]),
matrix + ".parilu" + output_suffix + "-u");
} else if (precond == "parilut") {
auto factory_parameter = gko::factorization::ParIlut<>::build();
if (argc >= 5) {
factory_parameter.with_iterations(
static_cast<gko::size_type>(std::stoi(argv[4])));
}
if (argc >= 6) {
factory_parameter.with_fill_in_limit(std::stod(argv[5]));
}
auto factory = factory_parameter.on(exec);
auto ilut = gko::as<gko::Composition<>>(
try_generate([&] { return factory->generate(mtx); }));
output(gko::as<gko::matrix::Csr<>>(ilut->get_operators()[0]),
matrix + ".parilut" + output_suffix + "-l");
output(gko::as<gko::matrix::Csr<>>(ilut->get_operators()[1]),
matrix + ".parilut" + output_suffix + "-u");
} else if (precond == "ilu-isai") {
auto fact_factory =
int sparsity_power = 1;
if (argc >= 5) {
sparsity_power = std::stoi(argv[4]);
}
auto factory =
.with_factorization(fact_factory)
.with_sparsity_power(sparsity_power))
.with_sparsity_power(sparsity_power))
.on(exec);
auto ilu_isai = try_generate([&] { return factory->generate(mtx); });
output(ilu_isai->get_l_solver()->get_approximate_inverse(),
matrix + ".ilu-isai" + output_suffix + "-l");
output(ilu_isai->get_u_solver()->get_approximate_inverse(),
matrix + ".ilu-isai" + output_suffix + "-u");
} else if (precond == "parilu-isai") {
auto fact_parameter = gko::factorization::ParIlu<>::build();
int sparsity_power = 1;
if (argc >= 5) {
fact_parameter.with_iterations(
static_cast<gko::size_type>(std::stoi(argv[4])));
}
if (argc >= 6) {
sparsity_power = std::stoi(argv[5]);
}
auto fact_factory = gko::share(fact_parameter.on(exec));
auto factory =
.with_factorization(fact_factory)
.with_sparsity_power(sparsity_power))
.with_sparsity_power(sparsity_power))
.on(exec);
auto ilu_isai = try_generate([&] { return factory->generate(mtx); });
output(ilu_isai->get_l_solver()->get_approximate_inverse(),
matrix + ".parilu-isai" + output_suffix + "-l");
output(ilu_isai->get_u_solver()->get_approximate_inverse(),
matrix + ".parilu-isai" + output_suffix + "-u");
} else if (precond == "parilut-isai") {
auto fact_parameter = gko::factorization::ParIlut<>::build();
int sparsity_power = 1;
if (argc >= 5) {
fact_parameter.with_iterations(
static_cast<gko::size_type>(std::stoi(argv[4])));
}
if (argc >= 6) {
fact_parameter.with_fill_in_limit(std::stod(argv[5]));
}
if (argc >= 7) {
sparsity_power = std::stoi(argv[6]);
}
auto fact_factory = gko::share(fact_parameter.on(exec));
auto factory =
.with_factorization(fact_factory)
.with_sparsity_power(sparsity_power))
.with_sparsity_power(sparsity_power))
.on(exec);
auto ilu_isai = try_generate([&] { return factory->generate(mtx); });
output(ilu_isai->get_l_solver()->get_approximate_inverse(),
matrix + ".parilut-isai" + output_suffix + "-l");
output(ilu_isai->get_u_solver()->get_approximate_inverse(),
matrix + ".parilut-isai" + output_suffix + "-u");
}
}
gko::Error::what
virtual const char * what() const noexcept override
Returns a human-readable string with a more detailed description of the error.
Definition: exception.hpp:75
gko::matrix::Csr
CSR is a matrix format which stores only the nonzero coefficients by compressing each row of the matr...
Definition: matrix.hpp:27
gko::max
constexpr T max(const T &x, const T &y)
Returns the larger of the arguments.
Definition: math.hpp:845
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:335
gko::log::profile_event_category::factory
LinOpFactory events.
gko::preconditioner::Isai
The Incomplete Sparse Approximate Inverse (ISAI) Preconditioner generates an approximate inverse matr...
Definition: isai.hpp:80
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
The Ginkgo namespace.
Definition: abstract_factory.hpp:20
gko::preconditioner::Jacobi
A block-Jacobi preconditioner is a block-diagonal linear operator, obtained by inverting the diagonal...
Definition: jacobi.hpp:187
gko::factorization::ParIlut
ParILUT is an incomplete threshold-based LU factorization which is computed in parallel.
Definition: par_ilut.hpp:73
gko::as
std::decay_t< T > * as(U *obj)
Performs polymorphic type conversion.
Definition: utils_helper.hpp:309
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:296
gko::ptr_param
This class is used for function parameters in the place of raw pointers.
Definition: utils_helper.hpp:43
gko::share
detail::shared_type< OwningPointer > share(OwningPointer &&p)
Marks the object pointed to by p as shared.
Definition: utils_helper.hpp:226
gko::factorization::ParIlu
ParILU is an incomplete LU factorization which is computed in parallel.
Definition: par_ilu.hpp:71
gko::WritableToMatrixData
A LinOp implementing this interface can write its data to a matrix_data structure.
Definition: lin_op.hpp:661
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:1345
gko::layout_type::coordinate
The matrix should be written as a sparse matrix in coordinate format.
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::Error
The Error class is used to report exceptional behaviour in library functions.
Definition: exception.hpp:58
gko::factorization::Ilu
Represents an incomplete LU factorization – ILU(0) – of a sparse matrix.
Definition: ilu.hpp:45
gko::preconditioner::Ilu
The Incomplete LU (ILU) preconditioner solves the equation for a given lower triangular matrix L,...
Definition: ilu.hpp:124