The custom matrix format example.
This example depends on simple-solver, poisson-solver, three-pt-stencil-solver, .
Introduction
This example solves a 1D Poisson equation:
![$ u : [0, 1] \rightarrow R\\ u'' = f\\ u(0) = u0\\ u(1) = u1 $](form_15.png)
using a finite difference method on an equidistant grid with K
discretization points (K
can be controlled with a command line parameter). The discretization is done via the second order Taylor polynomial:

For an equidistant grid with K "inner" discretization points
and step size
, the formula produces a system of linear equations

which is then solved using Ginkgo's implementation of the CG method preconditioned with block-Jacobi. It is also possible to specify on which executor Ginkgo will solve the system via the command line. The function
is set to
(making the solution
), but that can be changed in the main
function.
The intention of this example is to show how a custom linear operator can be created and integrated into Ginkgo to achieve performance benefits.
About the example
The commented program
#include <iostream>
#include <map>
#include <string>
#include <omp.h>
#include <ginkgo/ginkgo.hpp>
A CUDA kernel implementing the stencil, which will be used if running on the CUDA executor. Unfortunately, NVCC has serious problems interpreting some parts of Ginkgo's code, so the kernel has to be compiled separately.
extern void stencil_kernel(std::size_t size, const double *coefs,
const double *b, double *x);
A stencil matrix class representing the 3pt stencil linear operator. We include the gko::EnableLinOp mixin which implements the entire LinOp interface, except the two apply_impl methods, which get called inside the default implementation of apply (after argument verification) to perform the actual application of the linear operator. In addition, it includes the implementation of the entire PolymorphicObject interface.
It also includes the gko::EnableCreateMethod mixin which provides a default implementation of the static create method. This method will forward all its arguments to the constructor to create the object, and return an std::unique_ptr to the created object.
This constructor will be called by the create method. Here we initialize the coefficients of the stencil.
StencilMatrix(std::shared_ptr<const gko::Executor> exec,
double center = 2.0, double right = -1.0)
:
gko::EnableLinOp<StencilMatrix>(exec,
gko::dim<2>{size}),
coefficients(exec, {left, center, right})
{}
protected:
Here we implement the application of the linear operator, x = A * b. apply_impl will be called by the apply method, after the arguments have been moved to the correct executor and the operators checked for conforming sizes.
For simplicity, we assume that there is always only one right hand side and the stride of consecutive elements in the vectors is 1 (both of these are always true in this example).
we only implement the operator for dense RHS. gko::as will throw an exception if its argument is not Dense.
auto dense_b = gko::as<vec>(b);
auto dense_x = gko::as<vec>(x);
we need separate implementations depending on the executor, so we create an operation which maps the call to the correct implementation
stencil_operation(const coef_type &coefficients, const vec *b,
vec *x)
: coefficients{coefficients}, b{b}, x{x}
{}
OpenMP implementation
void run(std::shared_ptr<const gko::OmpExecutor>) const override
{
auto b_values = b->get_const_values();
auto x_values = x->get_values();
#pragma omp parallel for
for (std::size_t i = 0; i < x->
get_size()[0]; ++i) {
auto coefs = coefficients.get_const_data();
auto result = coefs[1] * b_values[i];
if (i > 0) {
result += coefs[0] * b_values[i - 1];
}
if (i < x->get_size()[0] - 1) {
result += coefs[2] * b_values[i + 1];
}
x_values[i] = result;
}
}
CUDA implementation
void run(std::shared_ptr<const gko::CudaExecutor>) const override
{
stencil_kernel(x->
get_size()[0], coefficients.get_const_data(),
b->get_const_values(), x->get_values());
}
We do not provide an implementation for reference executor. If not provided, Ginkgo will use the implementation for the OpenMP executor when calling it in the reference executor.
const coef_type &coefficients;
const vec *b;
vec *x;
};
this->get_executor()->run(
stencil_operation(coefficients, dense_b, dense_x));
}
There is also a version of the apply function which does the operation x = alpha * A * b + beta * x. This function is commonly used and can often be better optimized than implementing it using x = A * b. However, for simplicity, we will implement it exactly like that in this example.
{
auto dense_b = gko::as<vec>(b);
auto dense_x = gko::as<vec>(x);
auto tmp_x = dense_x->clone();
this->apply_impl(b,
lend(tmp_x));
dense_x->scale(beta);
dense_x->add_scaled(alpha,
lend(tmp_x));
}
private:
coef_type coefficients;
};
Creates a stencil matrix in CSR format for the given number of discretization points.
{
const auto discretization_points = matrix->
get_size()[0];
int pos = 0;
const double coefs[] = {-1, 2, -1};
row_ptrs[0] = pos;
for (int i = 0; i < discretization_points; ++i) {
for (auto ofs : {-1, 0, 1}) {
if (0 <= i + ofs && i + ofs < discretization_points) {
values[pos] = coefs[ofs + 1];
col_idxs[pos] = i + ofs;
++pos;
}
}
row_ptrs[i + 1] = pos;
}
}
Generates the RHS vector given f
and the boundary conditions.
template <typename Closure>
{
const auto discretization_points = rhs->
get_size()[0];
const auto h = 1.0 / (discretization_points + 1);
for (int i = 0; i < discretization_points; ++i) {
const auto xi = (i + 1) * h;
values[i] = -f(xi) * h * h;
}
values[0] += u0;
values[discretization_points - 1] += u1;
}
Prints the solution u
.
{
std::cout << u0 << '\n';
for (
int i = 0; i < u->
get_size()[0]; ++i) {
}
std::cout << u1 << std::endl;
}
Computes the 1-norm of the error given the computed u
and the correct solution function correct_u
.
template <typename Closure>
Closure correct_u)
{
const auto h = 1.0 / (discretization_points + 1);
auto error = 0.0;
for (int i = 0; i < discretization_points; ++i) {
using std::abs;
const auto xi = (i + 1) * h;
error +=
}
return error;
}
int main(int argc, char *argv[])
{
Some shortcuts
if (argc < 2) {
std::cerr << "Usage: " << argv[0] << " DISCRETIZATION_POINTS [executor]"
<< std::endl;
std::exit(-1);
}
Get number of discretization points
const unsigned int discretization_points =
argc >= 2 ? std::atoi(argv[1]) : 100u;
const auto executor_string = argc >= 3 ? argv[2] : "reference";
Figure out where to run the code
std::map<std::string, std::shared_ptr<gko::Executor>> exec_map{
{"omp", omp},
{"reference", gko::ReferenceExecutor::create()}};
executor where Ginkgo will perform the computation
const auto exec = exec_map.at(executor_string);
executor used by the application
const auto app_exec = exec_map["omp"];
problem:
auto correct_u = [](double x) { return x * x * x; };
auto f = [](double x) { return 6 * x; };
auto u0 = correct_u(0);
auto u1 = correct_u(1);
initialize vectors
auto rhs = vec::create(app_exec,
gko::dim<2>(discretization_points, 1));
generate_rhs(f, u0, u1,
lend(rhs));
auto u = vec::create(app_exec,
gko::dim<2>(discretization_points, 1));
for (
int i = 0; i < u->
get_size()[0]; ++i) {
}
Generate solver and solve the system
cg::build()
.with_criteria(gko::stop::Iteration::build()
.with_max_iters(discretization_points)
.on(exec),
.with_reduction_factor(1e-6)
.on(exec))
.on(exec)
notice how our custom StencilMatrix can be used in the same way as any built-in type
->generate(
StencilMatrix::create(exec, discretization_points, -1, 2, -1))
print_solution(u0, u1,
lend(u));
std::cout << "The average relative error is "
<< calculate_error(discretization_points,
lend(u), correct_u) /
discretization_points
<< std::endl;
}
Results
Comments about programming and debugging
The plain program
#include <iostream>
#include <map>
#include <string>
#include <omp.h>
#include <ginkgo/ginkgo.hpp>
extern void stencil_kernel(std::size_t size, const double *coefs,
const double *b, double *x);
public:
StencilMatrix(std::shared_ptr<const gko::Executor> exec,
double center = 2.0, double right = -1.0)
:
gko::EnableLinOp<StencilMatrix>(exec,
gko::dim<2>{size}),
coefficients(exec, {left, center, right})
{}
protected:
{
auto dense_b = gko::as<vec>(b);
auto dense_x = gko::as<vec>(x);
stencil_operation(const coef_type &coefficients, const vec *b,
vec *x)
: coefficients{coefficients}, b{b}, x{x}
{}
void run(std::shared_ptr<const gko::OmpExecutor>) const override
{
auto b_values = b->get_const_values();
auto x_values = x->get_values();
#pragma omp parallel for
for (std::size_t i = 0; i < x->
get_size()[0]; ++i) {
auto coefs = coefficients.get_const_data();
auto result = coefs[1] * b_values[i];
if (i > 0) {
result += coefs[0] * b_values[i - 1];
}
if (i < x->get_size()[0] - 1) {
result += coefs[2] * b_values[i + 1];
}
x_values[i] = result;
}
}
void run(std::shared_ptr<const gko::CudaExecutor>) const override
{
stencil_kernel(x->
get_size()[0], coefficients.get_const_data(),
b->get_const_values(), x->get_values());
}
const coef_type &coefficients;
const vec *b;
vec *x;
};
this->get_executor()->run(
stencil_operation(coefficients, dense_b, dense_x));
}
{
auto dense_b = gko::as<vec>(b);
auto dense_x = gko::as<vec>(x);
auto tmp_x = dense_x->clone();
this->apply_impl(b,
lend(tmp_x));
dense_x->scale(beta);
dense_x->add_scaled(alpha,
lend(tmp_x));
}
private:
coef_type coefficients;
};
{
const auto discretization_points = matrix->
get_size()[0];
int pos = 0;
const double coefs[] = {-1, 2, -1};
row_ptrs[0] = pos;
for (int i = 0; i < discretization_points; ++i) {
for (auto ofs : {-1, 0, 1}) {
if (0 <= i + ofs && i + ofs < discretization_points) {
values[pos] = coefs[ofs + 1];
col_idxs[pos] = i + ofs;
++pos;
}
}
row_ptrs[i + 1] = pos;
}
}
template <typename Closure>
{
const auto discretization_points = rhs->
get_size()[0];
const auto h = 1.0 / (discretization_points + 1);
for (int i = 0; i < discretization_points; ++i) {
const auto xi = (i + 1) * h;
values[i] = -f(xi) * h * h;
}
values[0] += u0;
values[discretization_points - 1] += u1;
}
{
std::cout << u0 << '\n';
for (
int i = 0; i < u->
get_size()[0]; ++i) {
}
std::cout << u1 << std::endl;
}
template <typename Closure>
Closure correct_u)
{
const auto h = 1.0 / (discretization_points + 1);
auto error = 0.0;
for (int i = 0; i < discretization_points; ++i) {
using std::abs;
const auto xi = (i + 1) * h;
error +=
}
return error;
}
int main(int argc, char *argv[])
{
if (argc < 2) {
std::cerr << "Usage: " << argv[0] << " DISCRETIZATION_POINTS [executor]"
<< std::endl;
std::exit(-1);
}
const unsigned int discretization_points =
argc >= 2 ? std::atoi(argv[1]) : 100u;
const auto executor_string = argc >= 3 ? argv[2] : "reference";
std::map<std::string, std::shared_ptr<gko::Executor>> exec_map{
{"omp", omp},
{"reference", gko::ReferenceExecutor::create()}};
const auto exec = exec_map.at(executor_string);
const auto app_exec = exec_map["omp"];
auto correct_u = [](double x) { return x * x * x; };
auto f = [](double x) { return 6 * x; };
auto u0 = correct_u(0);
auto u1 = correct_u(1);
auto rhs = vec::create(app_exec,
gko::dim<2>(discretization_points, 1));
generate_rhs(f, u0, u1,
lend(rhs));
auto u = vec::create(app_exec,
gko::dim<2>(discretization_points, 1));
for (
int i = 0; i < u->
get_size()[0]; ++i) {
}
cg::build()
.with_criteria(gko::stop::Iteration::build()
.with_max_iters(discretization_points)
.on(exec),
.with_reduction_factor(1e-6)
.on(exec))
.on(exec)
->generate(
StencilMatrix::create(exec, discretization_points, -1, 2, -1))
print_solution(u0, u1,
lend(u));
std::cout << "The average relative error is "
<< calculate_error(discretization_points,
lend(u), correct_u) /
discretization_points
<< std::endl;
}