Ginkgo  Generated from pipelines/130473384 branch based on develop. Ginkgo version 1.1.1
A numerical linear algebra library targeting many-core architectures
The poisson-solver program

The poisson solver example..

This example depends on simple-solver.

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

Introduction

This example solves a 1D Poisson equation:

$ u : [0, 1] \rightarrow R\\ u'' = f\\ u(0) = u0\\ u(1) = u1 $

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:

$ u(x + h) = u(x) - u'(x)h + 1/2 u''(x)h^2 + O(h^3)\\ u(x - h) = u(x) + u'(x)h + 1/2 u''(x)h^2 + O(h^3) / +\\ ---------------------- \\ -u(x - h) + 2u(x) + -u(x + h) = -f(x)h^2 + O(h^3) $

For an equidistant grid with K "inner" discretization points $x1, ..., xk, $and step size $ h = 1 / (K + 1)$, the formula produces a system of linear equations

$ 2u_1 - u_2 = -f_1 h^2 + u0\\ -u_(k-1) + 2u_k - u_(k+1) = -f_k h^2, k = 2, ..., K - 1\\ -u_(K-1) + 2u_K = -f_K h^2 + u1\\ $

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 $`f` $is set to $`f(x) = 6x`$ (making the solution $`u(x) = x^3`$), but that can be changed in the main function.

The intention of the example is to show how Ginkgo can be used to build an application solving a real-world problem, which includes a solution of a large, sparse linear system as a component.

About the example

The commented program

#include <ginkgo/ginkgo.hpp>
#include <iostream>
#include <map>
#include <string>
#include <vector>

Creates a stencil matrix in CSR format for the given number of discretization points.

void generate_stencil_matrix(gko::matrix::Csr<> *matrix)
{
const auto discretization_points = matrix->get_size()[0];
auto row_ptrs = matrix->get_row_ptrs();
auto col_idxs = matrix->get_col_idxs();
auto values = matrix->get_values();
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>
void generate_rhs(Closure f, double u0, double u1, gko::matrix::Dense<> *rhs)
{
const auto discretization_points = rhs->get_size()[0];
auto values = rhs->get_values();
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.

void print_solution(double u0, double u1, const gko::matrix::Dense<> *u)
{
std::cout << u0 << '\n';
for (int i = 0; i < u->get_size()[0]; ++i) {
std::cout << u->get_const_values()[i] << '\n';
}
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>
double calculate_error(int discretization_points, const gko::matrix::Dense<> *u,
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 +=
abs(u->get_const_values()[i] - correct_u(xi)) / abs(correct_u(xi));
}
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]) : 100;
const auto executor_string = argc >= 3 ? argv[2] : "reference";

Figure out where to run the code

const auto omp = gko::OmpExecutor::create();
std::map<std::string, std::shared_ptr<gko::Executor>> exec_map{
{"omp", omp},
{"cuda", gko::CudaExecutor::create(0, omp)},
{"hip", gko::HipExecutor::create(0, omp)},
{"reference", gko::ReferenceExecutor::create()}};

executor where Ginkgo will perform the computation

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

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 matrix and vectors

auto matrix = mtx::create(app_exec, gko::dim<2>(discretization_points),
3 * discretization_points - 2);
generate_stencil_matrix(lend(matrix));
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) {
u->get_values()[i] = 0.0;
}

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))
.with_preconditioner(bj::build().on(exec))
.on(exec)
->generate(clone(exec, matrix)) // copy the matrix to the executor
->apply(lend(rhs), lend(u));
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

This is the expected output:

0
0.00010798
0.000863838
0.00291545
0.0069107
0.0134975
0.0233236
0.037037
0.0552856
0.0787172
0.10798
0.143721
0.186589
0.237231
0.296296
0.364431
0.442285
0.530504
0.629738
0.740633
0.863838
1
The average relative error is 1.87318e-15

Comments about programming and debugging

The plain program

/*******************************<GINKGO LICENSE>******************************
Copyright (c) 2017-2020, 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 <iostream>
#include <map>
#include <string>
#include <vector>
void generate_stencil_matrix(gko::matrix::Csr<> *matrix)
{
const auto discretization_points = matrix->get_size()[0];
auto row_ptrs = matrix->get_row_ptrs();
auto col_idxs = matrix->get_col_idxs();
auto values = matrix->get_values();
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>
void generate_rhs(Closure f, double u0, double u1, gko::matrix::Dense<> *rhs)
{
const auto discretization_points = rhs->get_size()[0];
auto values = rhs->get_values();
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;
}
void print_solution(double u0, double u1, const gko::matrix::Dense<> *u)
{
std::cout << u0 << '\n';
for (int i = 0; i < u->get_size()[0]; ++i) {
std::cout << u->get_const_values()[i] << '\n';
}
std::cout << u1 << std::endl;
}
template <typename Closure>
double calculate_error(int discretization_points, const gko::matrix::Dense<> *u,
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 +=
abs(u->get_const_values()[i] - correct_u(xi)) / abs(correct_u(xi));
}
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]) : 100;
const auto executor_string = argc >= 3 ? argv[2] : "reference";
const auto omp = gko::OmpExecutor::create();
std::map<std::string, std::shared_ptr<gko::Executor>> exec_map{
{"omp", omp},
{"cuda", gko::CudaExecutor::create(0, omp)},
{"hip", gko::HipExecutor::create(0, omp)},
{"reference", gko::ReferenceExecutor::create()}};
const auto exec = exec_map.at(executor_string); // throws if not valid
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 matrix = mtx::create(app_exec, gko::dim<2>(discretization_points),
3 * discretization_points - 2);
generate_stencil_matrix(lend(matrix));
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) {
u->get_values()[i] = 0.0;
}
cg::build()
.with_criteria(gko::stop::Iteration::build()
.with_max_iters(discretization_points)
.on(exec),
.with_reduction_factor(1e-6)
.on(exec))
.with_preconditioner(bj::build().on(exec))
.on(exec)
->generate(clone(exec, matrix)) // copy the matrix to the executor
->apply(lend(rhs), lend(u));
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;
}
gko::abs
constexpr T abs(const T &x)
Returns the absolute value of the object.
Definition: math.hpp:572
gko::matrix::Csr::get_col_idxs
index_type * get_col_idxs() noexcept
Returns the column indexes of the matrix.
Definition: csr.hpp:647
gko::matrix::Csr
CSR is a matrix format which stores only the nonzero coefficients by compressing each row of the matr...
Definition: coo.hpp:51
gko::HipExecutor::create
static std::shared_ptr< HipExecutor > create(int device_id, std::shared_ptr< Executor > master)
Creates a new HipExecutor.
gko::matrix::Dense
Dense is a matrix format which explicitly stores all values of the matrix.
Definition: coo.hpp:55
gko::CudaExecutor::create
static std::shared_ptr< CudaExecutor > create(int device_id, std::shared_ptr< Executor > master)
Creates a new CudaExecutor.
gko::OmpExecutor::create
static std::shared_ptr< OmpExecutor > create()
Creates a new OmpExecutor.
Definition: executor.hpp:775
gko::stop::ResidualNormReduction
The ResidualNormReduction class is a stopping criterion which stops the iteration process when the re...
Definition: residual_norm_reduction.hpp:64
gko::clone
detail::cloned_type< Pointer > clone(const Pointer &p)
Creates a unique clone of the object pointed to by p.
Definition: utils.hpp:160
gko::matrix::Dense::get_values
value_type * get_values() noexcept
Returns a pointer to the array of values of the matrix.
Definition: dense.hpp:239
gko::preconditioner::Jacobi
A block-Jacobi preconditioner is a block-diagonal linear operator, obtained by inverting the diagonal...
Definition: jacobi.hpp:207
gko::matrix::Csr::get_row_ptrs
index_type * get_row_ptrs() noexcept
Returns the row pointers of the matrix.
Definition: csr.hpp:666
gko::lend
std::enable_if< detail::have_ownership_s< Pointer >::value, detail::pointee< Pointer > * >::type lend(const Pointer &p)
Returns a non-owning (plain) pointer to the object pointed to by p.
Definition: utils.hpp:253
gko::dim< 2 >
gko::solver::Cg
CG or the conjugate gradient method is an iterative type Krylov subspace method which is suitable for...
Definition: cg.hpp:72
gko::matrix::Dense::get_const_values
const value_type * get_const_values() const noexcept
Returns a pointer to the array of values of the matrix.
Definition: dense.hpp:248
gko::LinOp::get_size
const dim< 2 > & get_size() const noexcept
Returns the size of the operator.
Definition: lin_op.hpp:221
gko::matrix::Csr::get_values
value_type * get_values() noexcept
Returns the values of the matrix.
Definition: csr.hpp:628