Ginkgo  Generated from pipelines/1570051012 branch based on develop. Ginkgo version 1.9.0
A numerical linear algebra library targeting many-core architectures
The file-config-solver program

The file config solver example..

This example depends on simple-solver.

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

This example shows how to use file to configure solver.

In this example, we first read in a matrix from a file. Then we read the config file to configure the solver. The example application reports the generation and run time of the solver.

The commented program

const std::string matrix_path = argc >= 4 ? argv[3] : "data/A.mtx";

Figure out where to run the code

std::map<std::string, std::function<std::shared_ptr<gko::Executor>()>>
exec_map{
{"omp", [] { return gko::OmpExecutor::create(); }},
{"cuda",
[] {
}},
{"hip",
[] {
}},
{"dpcpp",
[] {
0, gko::ReferenceExecutor::create());
}},
{"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 = share(gko::read<mtx>(std::ifstream(matrix_path), exec));

Create RHS as 1 and initial guess as 0

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

Calculate initial residual by overwriting b

auto one = gko::initialize<vec>({1.0}, exec);
auto neg_one = gko::initialize<vec>({-1.0}, exec);
auto initres = gko::initialize<vec>({0.0}, exec);
A->apply(one, x, neg_one, b);
b->compute_norm2(initres);

Copy b again

b->copy_from(host_b);

Read the json config file to configure the ginkgo solver. The following files, which are mapped to corresponding examples, are available cg.json: simple-solver blockjacobi-cg.json: preconditioned-solver ir.json: iterative-refinement parilu.json: ilu-preconditioned-solver (by using factoization parameter directly) pgm-multigrid-cg.json: multigrid-preconditioned-solver (set min_coarse_rows additionally due to this small example matrix) mixed-pgm-multigrid-cg.json: mixed-multigrid-preconditioned-solver (assuming there are always more than one level)

auto config = gko::ext::config::parse_json_file(configfile);

Create the registry, which allows passing the existing data into config This example does not use existing data.

auto reg = gko::config::registry();

Create the default type descriptor, which gives the default common type (value/index) for solver generation. If the solver does not specify value type, the solver will use these types.

auto td = gko::config::make_type_descriptor<ValueType, IndexType>();

generate the linopfactory on the given executors

auto solver_gen = gko::config::parse(config, reg, td).on(exec);

Create solver

const auto gen_tic = std::chrono::steady_clock::now();
auto solver = solver_gen->generate(A);
exec->synchronize();
const auto gen_toc = std::chrono::steady_clock::now();
const auto gen_time = std::chrono::duration_cast<std::chrono::milliseconds>(
gen_toc - gen_tic);

Add logger

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

Solve system

exec->synchronize();
const auto tic = std::chrono::steady_clock::now();
solver->apply(b, x);
exec->synchronize();
const auto toc = std::chrono::steady_clock::now();
const auto time =
std::chrono::duration_cast<std::chrono::milliseconds>(toc - tic);

Print out the solver config

std::cout << "Config file: " << configfile << std::endl;
std::ifstream f(configfile);
std::cout << f.rdbuf() << std::endl;

Calculate residual

auto res = gko::as<vec>(logger->get_residual_norm());
std::cout << "Initial residual norm sqrt(r^T r): \n";
write(std::cout, initres);
std::cout << "Final residual norm sqrt(r^T r): \n";
write(std::cout, res);

Print solver statistics

std::cout << "Solver iteration count: " << logger->get_num_iterations()
<< std::endl;
std::cout << "Solver generation time [ms]: "
<< static_cast<double>(gen_time.count()) << std::endl;
std::cout << "Solver execution time [ms]: "
<< static_cast<double>(time.count()) << std::endl;
std::cout << "Solver execution time per iteration[ms]: "
<< static_cast<double>(time.count()) /
logger->get_num_iterations()
<< std::endl;
}

Results

This is the expected output:

Config file: config/cg.json
{
"type": "solver::Cg",
"criteria": [
{
"type": "Iteration",
"max_iters": 20
},
{
"type": "ResidualNorm",
"reduction_factor": 1e-7
}
]
}
Initial residual norm sqrt(r^T r):
%%MatrixMarket matrix array real general
1 1
25.9808
Final residual norm sqrt(r^T r):
%%MatrixMarket matrix array real general
1 1
58.79
Solver iteration count: 20
Solver generation time [ms]: 0.065244
Solver execution time [ms]: 0.793764
Solver execution time per iteration[ms]: 0.0396882

Comments about programming and debugging

The plain program

#include <chrono>
#include <fstream>
#include <iomanip>
#include <iostream>
#include <map>
#include <string>
#include <ginkgo/ginkgo.hpp>
#include <ginkgo/extensions/config/json_config.hpp>
int main(int argc, char* argv[])
{
using ValueType = double;
using IndexType = int;
std::cout << gko::version_info::get() << std::endl;
std::cout << argv[0] << " executor configfile matrix" << std::endl;
const auto executor_string = argc >= 2 ? argv[1] : "reference";
const auto configfile = argc >= 3 ? argv[2] : "config/cg.json";
const std::string matrix_path = argc >= 4 ? argv[3] : "data/A.mtx";
std::map<std::string, std::function<std::shared_ptr<gko::Executor>()>>
exec_map{
{"omp", [] { return gko::OmpExecutor::create(); }},
{"cuda",
[] {
}},
{"hip",
[] {
}},
{"dpcpp",
[] {
0, gko::ReferenceExecutor::create());
}},
{"reference", [] { return gko::ReferenceExecutor::create(); }}};
const auto exec = exec_map.at(executor_string)(); // throws if not valid
auto A = share(gko::read<mtx>(std::ifstream(matrix_path), exec));
gko::size_type size = A->get_size()[0];
auto host_x = vec::create(exec->get_master(), gko::dim<2>(size, 1));
auto host_b = vec::create(exec->get_master(), gko::dim<2>(size, 1));
for (auto i = 0; i < size; i++) {
host_x->at(i, 0) = 0.;
host_b->at(i, 0) = 1.;
}
auto x = vec::create(exec);
auto b = vec::create(exec);
x->copy_from(host_x);
b->copy_from(host_b);
auto one = gko::initialize<vec>({1.0}, exec);
auto neg_one = gko::initialize<vec>({-1.0}, exec);
auto initres = gko::initialize<vec>({0.0}, exec);
A->apply(one, x, neg_one, b);
b->compute_norm2(initres);
b->copy_from(host_b);
auto config = gko::ext::config::parse_json_file(configfile);
auto reg = gko::config::registry();
auto td = gko::config::make_type_descriptor<ValueType, IndexType>();
auto solver_gen = gko::config::parse(config, reg, td).on(exec);
const auto gen_tic = std::chrono::steady_clock::now();
auto solver = solver_gen->generate(A);
exec->synchronize();
const auto gen_toc = std::chrono::steady_clock::now();
const auto gen_time = std::chrono::duration_cast<std::chrono::milliseconds>(
gen_toc - gen_tic);
std::shared_ptr<const gko::log::Convergence<ValueType>> logger =
solver->add_logger(logger);
exec->synchronize();
const auto tic = std::chrono::steady_clock::now();
solver->apply(b, x);
exec->synchronize();
const auto toc = std::chrono::steady_clock::now();
const auto time =
std::chrono::duration_cast<std::chrono::milliseconds>(toc - tic);
std::cout << "Config file: " << configfile << std::endl;
std::ifstream f(configfile);
std::cout << f.rdbuf() << std::endl;
auto res = gko::as<vec>(logger->get_residual_norm());
std::cout << "Initial residual norm sqrt(r^T r): \n";
write(std::cout, initres);
std::cout << "Final residual norm sqrt(r^T r): \n";
write(std::cout, res);
std::cout << "Solver iteration count: " << logger->get_num_iterations()
<< std::endl;
std::cout << "Solver generation time [ms]: "
<< static_cast<double>(gen_time.count()) << std::endl;
std::cout << "Solver execution time [ms]: "
<< static_cast<double>(time.count()) << std::endl;
std::cout << "Solver execution time per iteration[ms]: "
<< static_cast<double>(time.count()) /
logger->get_num_iterations()
<< std::endl;
}
gko::matrix::Csr
CSR is a matrix format which stores only the nonzero coefficients by compressing each row of the matr...
Definition: matrix.hpp:28
gko::log::profile_event_category::solver
Solver events.
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:73
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:19
gko::size_type
std::size_t size_type
Integral type used for allocation quantities.
Definition: types.hpp:89
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:139
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:295
gko::share
detail::shared_type< OwningPointer > share(OwningPointer &&p)
Marks the object pointed to by p as shared.
Definition: utils_helper.hpp:224
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:1396
gko::config::registry
This class stores additional context for creating Ginkgo objects from configuration files.
Definition: registry.hpp:167
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:871
gko::one
constexpr T one()
Returns the multiplicative identity for T.
Definition: math.hpp:632