Ginkgo Generated from branch based on main. Ginkgo version 1.9.0
A numerical linear algebra library targeting many-core architectures
 
Loading...
Searching...
No Matches
The three-pt-stencil-solver program

The 3-point stencil example.

This example depends on simple-solver, poisson-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 integrated into existing software - the generate_stencil_matrix, generate_rhs, print_solution, compute_error and main function do not reference Ginkgo at all (i.e. they could have been there before the application developer decided to use Ginkgo, and the only part where Ginkgo is introduced is inside the solve_system function.

About the example

The commented program

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 integrated into
existing software - the `generate_stencil_matrix`, `generate_rhs`,
`print_solution`, `compute_error` and `main` function do not reference Ginkgo at
all (i.e. they could have been there before the application developer decided to
use Ginkgo, and the only part where Ginkgo is introduced is inside the
`solve_system` function.
*****************************<DESCRIPTION>********************************** /
#include <iostream>
#include <map>
#include <string>
#include <vector>
#include <ginkgo/ginkgo.hpp>

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

template <typename ValueType, typename IndexType>
void generate_stencil_matrix(IndexType discretization_points,
IndexType* row_ptrs, IndexType* col_idxs,
ValueType* values)
{
IndexType pos = 0;
const ValueType coefs[] = {-1, 2, -1};
row_ptrs[0] = pos;
for (IndexType 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, typename ValueType, typename IndexType>
void generate_rhs(IndexType discretization_points, Closure f, ValueType u0,
ValueType u1, ValueType* rhs)
{
const ValueType h = 1.0 / (discretization_points + 1);
for (IndexType i = 0; i < discretization_points; ++i) {
const ValueType xi = ValueType(i + 1) * h;
rhs[i] = -f(xi) * h * h;
}
rhs[0] += u0;
rhs[discretization_points - 1] += u1;
}

Prints the solution u.

template <typename ValueType, typename IndexType>
void print_solution(IndexType discretization_points, ValueType u0, ValueType u1,
const ValueType* u)
{
std::cout << u0 << '\n';
for (IndexType i = 0; i < discretization_points; ++i) {
std::cout << u[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, typename ValueType, typename IndexType>
gko::remove_complex<ValueType> calculate_error(IndexType discretization_points,
const ValueType* u,
Closure correct_u)
{
const ValueType h = 1.0 / (discretization_points + 1);
for (IndexType i = 0; i < discretization_points; ++i) {
using std::abs;
const ValueType xi = ValueType(i + 1) * h;
error += abs(u[i] - correct_u(xi)) / abs(correct_u(xi));
}
return error;
}
template <typename ValueType, typename IndexType>
void solve_system(const std::string& executor_string,
IndexType discretization_points, IndexType* row_ptrs,
IndexType* col_idxs, ValueType* values, ValueType* rhs,
ValueType* u, gko::remove_complex<ValueType> reduction_factor)
{
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:260

Some shortcuts

using val_array = gko::array<ValueType>;
using idx_array = gko::array<IndexType>;
const auto& dp = discretization_points;
An array is a container which encapsulates fixed-sized arrays, stored on the Executor tied to the arr...
Definition array.hpp:166
CSR is a matrix format which stores only the nonzero coefficients by compressing each row of the matr...
Definition csr.hpp:121
Dense is a matrix format which explicitly stores all values of the matrix.
Definition dense.hpp:117
A block-Jacobi preconditioner is a block-diagonal linear operator, obtained by inverting the diagonal...
Definition jacobi.hpp:190
CG or the conjugate gradient method is an iterative type Krylov subspace method which is suitable for...
Definition cg.hpp:50

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",
[] {
}},
{"reference", [] { return gko::ReferenceExecutor::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.
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.
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.
static std::shared_ptr< OmpExecutor > create(std::shared_ptr< CpuAllocatorBase > alloc=std::make_shared< CpuAllocator >())
Creates a new OmpExecutor.
Definition executor.hpp:1396

executor where Ginkgo will perform the computation

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

executor where the application initialized the data

const auto app_exec = exec->get_master();

Tell Ginkgo to use the data in our application

Matrix: we have to set the executor of the matrix to the one where we want SpMVs to run (in this case exec). When creating array views, we have to specify the executor where the data is (in this case app_exec).

If the two do not match, Ginkgo will automatically create a copy of the data on exec (however, it will not copy the data back once it is done

  • here this is not important since we are not modifying the matrix).
auto matrix = mtx::create(exec, gko::dim<2>(dp),
val_array::view(app_exec, 3 * dp - 2, values),
idx_array::view(app_exec, 3 * dp - 2, col_idxs),
idx_array::view(app_exec, dp + 1, row_ptrs));
A type representing the dimensions of a multidimensional object.
Definition dim.hpp:26

RHS: similar to matrix

auto b = vec::create(exec, gko::dim<2>(dp, 1),
val_array::view(app_exec, dp, rhs), 1);

Solution: we have to be careful here - if the executors are different, once we compute the solution the array will not be automatically copied back to the original memory locations. Fortunately, whenever apply is called on a linear operator (e.g. matrix, solver) the arguments automatically get copied to the executor where the operator is, and copied back once the operation is completed. Thus, in this case, we can just define the solution on app_exec, and it will be automatically transferred to/from exec if needed.

auto x = vec::create(app_exec, gko::dim<2>(dp, 1),
val_array::view(app_exec, dp, u), 1);

Generate solver

auto solver_gen =
cg::build()
.with_criteria(gko::stop::Iteration::build().with_max_iters(
gko::stop::ResidualNorm<ValueType>::build()
.with_reduction_factor(reduction_factor))
.with_preconditioner(bj::build())
.on(exec);
auto solver = solver_gen->generate(gko::give(matrix));
std::size_t size_type
Integral type used for allocation quantities.
Definition types.hpp:89
std::remove_reference< OwningPointer >::type && give(OwningPointer &&p)
Marks that the object pointed to by p can be given to the callee.
Definition utils_helper.hpp:247

Solve system

solver->apply(b, x);
}
int main(int argc, char* argv[])
{
using ValueType = double;
using IndexType = int;

Print version information

std::cout << gko::version_info::get() << std::endl;
if (argc == 2 && std::string(argv[1]) == "--help") {
std::cerr << "Usage: " << argv[0]
<< " [executor] [DISCRETIZATION_POINTS]" << std::endl;
std::exit(-1);
}
const auto executor_string = argc >= 2 ? argv[1] : "reference";
const IndexType discretization_points =
argc >= 3 ? std::atoi(argv[2]) : 100;
static const version_info & get()
Returns an instance of version_info.
Definition version.hpp:139

problem:

auto correct_u = [](ValueType x) { return x * x * x; };
auto f = [](ValueType x) { return ValueType(6) * x; };
auto u0 = correct_u(0);
auto u1 = correct_u(1);

matrix

std::vector<IndexType> row_ptrs(discretization_points + 1);
std::vector<IndexType> col_idxs(3 * discretization_points - 2);
std::vector<ValueType> values(3 * discretization_points - 2);

right hand side

std::vector<ValueType> rhs(discretization_points);

solution

std::vector<ValueType> u(discretization_points, 0.0);
const gko::remove_complex<ValueType> reduction_factor = 1e-7;
generate_stencil_matrix(discretization_points, row_ptrs.data(),
col_idxs.data(), values.data());

looking for solution u = x^3: f = 6x, u(0) = 0, u(1) = 1

generate_rhs(discretization_points, f, u0, u1, rhs.data());
solve_system(executor_string, discretization_points, row_ptrs.data(),
col_idxs.data(), values.data(), rhs.data(), u.data(),
reduction_factor);

Uncomment to print the solution print_solution<ValueType, IndexType>(discretization_points, 0, 1, u.data());

std::cout << "The average relative error is "
<< calculate_error(discretization_points, u.data(), correct_u) /
discretization_points
<< std::endl;
}

Results

This is the expected output:

The average relative error is 2.52236e-11

Comments about programming and debugging

The plain program

/*****************************<DESCRIPTION>***********************************
This example solves a 1D Poisson equation:
u : [0, 1] -> 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 integrated into
existing software - the `generate_stencil_matrix`, `generate_rhs`,
`print_solution`, `compute_error` and `main` function do not reference Ginkgo at
all (i.e. they could have been there before the application developer decided to
use Ginkgo, and the only part where Ginkgo is introduced is inside the
`solve_system` function.
*****************************<DESCRIPTION>**********************************/
#include <iostream>
#include <map>
#include <string>
#include <vector>
#include <ginkgo/ginkgo.hpp>
template <typename ValueType, typename IndexType>
void generate_stencil_matrix(IndexType discretization_points,
IndexType* row_ptrs, IndexType* col_idxs,
ValueType* values)
{
IndexType pos = 0;
const ValueType coefs[] = {-1, 2, -1};
row_ptrs[0] = pos;
for (IndexType 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, typename ValueType, typename IndexType>
void generate_rhs(IndexType discretization_points, Closure f, ValueType u0,
ValueType u1, ValueType* rhs)
{
const ValueType h = 1.0 / (discretization_points + 1);
for (IndexType i = 0; i < discretization_points; ++i) {
const ValueType xi = ValueType(i + 1) * h;
rhs[i] = -f(xi) * h * h;
}
rhs[0] += u0;
rhs[discretization_points - 1] += u1;
}
template <typename ValueType, typename IndexType>
void print_solution(IndexType discretization_points, ValueType u0, ValueType u1,
const ValueType* u)
{
std::cout << u0 << '\n';
for (IndexType i = 0; i < discretization_points; ++i) {
std::cout << u[i] << '\n';
}
std::cout << u1 << std::endl;
}
template <typename Closure, typename ValueType, typename IndexType>
gko::remove_complex<ValueType> calculate_error(IndexType discretization_points,
const ValueType* u,
Closure correct_u)
{
const ValueType h = 1.0 / (discretization_points + 1);
for (IndexType i = 0; i < discretization_points; ++i) {
using std::abs;
const ValueType xi = ValueType(i + 1) * h;
error += abs(u[i] - correct_u(xi)) / abs(correct_u(xi));
}
return error;
}
template <typename ValueType, typename IndexType>
void solve_system(const std::string& executor_string,
IndexType discretization_points, IndexType* row_ptrs,
IndexType* col_idxs, ValueType* values, ValueType* rhs,
ValueType* u, gko::remove_complex<ValueType> reduction_factor)
{
using val_array = gko::array<ValueType>;
using idx_array = gko::array<IndexType>;
const auto& dp = discretization_points;
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
const auto app_exec = exec->get_master();
auto matrix = mtx::create(exec, gko::dim<2>(dp),
val_array::view(app_exec, 3 * dp - 2, values),
idx_array::view(app_exec, 3 * dp - 2, col_idxs),
idx_array::view(app_exec, dp + 1, row_ptrs));
auto b = vec::create(exec, gko::dim<2>(dp, 1),
val_array::view(app_exec, dp, rhs), 1);
auto x = vec::create(app_exec, gko::dim<2>(dp, 1),
val_array::view(app_exec, dp, u), 1);
auto solver_gen =
cg::build()
.with_criteria(gko::stop::Iteration::build().with_max_iters(
gko::stop::ResidualNorm<ValueType>::build()
.with_reduction_factor(reduction_factor))
.with_preconditioner(bj::build())
.on(exec);
auto solver = solver_gen->generate(gko::give(matrix));
solver->apply(b, x);
}
int main(int argc, char* argv[])
{
using ValueType = double;
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] [DISCRETIZATION_POINTS]" << std::endl;
std::exit(-1);
}
const auto executor_string = argc >= 2 ? argv[1] : "reference";
const IndexType discretization_points =
argc >= 3 ? std::atoi(argv[2]) : 100;
auto correct_u = [](ValueType x) { return x * x * x; };
auto f = [](ValueType x) { return ValueType(6) * x; };
auto u0 = correct_u(0);
auto u1 = correct_u(1);
std::vector<IndexType> row_ptrs(discretization_points + 1);
std::vector<IndexType> col_idxs(3 * discretization_points - 2);
std::vector<ValueType> values(3 * discretization_points - 2);
std::vector<ValueType> rhs(discretization_points);
std::vector<ValueType> u(discretization_points, 0.0);
const gko::remove_complex<ValueType> reduction_factor = 1e-7;
generate_stencil_matrix(discretization_points, row_ptrs.data(),
col_idxs.data(), values.data());
generate_rhs(discretization_points, f, u0, u1, rhs.data());
solve_system(executor_string, discretization_points, row_ptrs.data(),
col_idxs.data(), values.data(), rhs.data(), u.data(),
reduction_factor);
std::cout << "The average relative error is "
<< calculate_error(discretization_points, u.data(), correct_u) /
discretization_points
<< std::endl;
}
@ solver
Solver events.
Definition profiler_hook.hpp:34
@ rhs
the input is right hand side
Definition solver_base.hpp:41
constexpr std::enable_if_t<!is_complex_s< T >::value, T > abs(const T &x)
Returns the absolute value of the object.
Definition math.hpp:931