Ginkgo  Generated from pipelines/1068515030 branch based on master. Ginkgo version 1.7.0
A numerical linear algebra library targeting many-core architectures
The kokkos_assembly program

The Kokkos assembly example..

This example depends on simple-solver, poisson-solver, three-pt-stencil-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 resulting CSR matrix is assembled using Kokkos kernels. This example show how to use Ginkgo data with Kokkos kernels.

Notes

If this example is built as part of Ginkgo, it is advised to configure Ginkgo with -DGINKGO_WITH_CCACHE=OFF to prevent incompabilities with Kokkos' compiler wrapper for nvcc.

The commented program

#include <iostream>
#include <map>
#include <string>
#include <omp.h>
#include <Kokkos_Core.hpp>
#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(gko::matrix::Csr<ValueType, IndexType>* matrix)
{
auto exec = matrix->get_executor();
const auto discretization_points = matrix->get_size()[0];

Over-allocate storage for the matrix elements. Each row has 3 entries, except for the first and last one. To handle each row uniformly, we allocate space for 3x the number of rows.

discretization_points * 3);

Create Kokkos views on Ginkgo data.

Kokkos::View<IndexType*> v_row_idxs(md.get_row_idxs(), md.get_num_elems());
Kokkos::View<IndexType*> v_col_idxs(md.get_col_idxs(), md.get_num_elems());
Kokkos::View<ValueType*> v_values(md.get_values(), md.get_num_elems());

Create the matrix entries. This also creates zero entries for the first and second row to handle all rows uniformly.

Kokkos::parallel_for(
"generate_stencil_matrix", md.get_num_elems(), KOKKOS_LAMBDA(int i) {
const ValueType coefs[] = {-1, 2, -1};
auto ofs = static_cast<IndexType>((i % 3) - 1);
auto row = static_cast<IndexType>(i / 3);
auto col = row + ofs;

To prevent branching, a mask is used to set the entry to zero, if the column is out-of-bounds

auto mask =
static_cast<IndexType>(0 <= col && col < discretization_points);
v_row_idxs[i] = mask * row;
v_col_idxs[i] = mask * col;
v_values[i] = mask * coefs[ofs + 1];
});

Add up duplicate (zero) entries.

md.sum_duplicates();

Build Csr matrix.

matrix->read(std::move(md));
}

Generates the RHS vector given f and the boundary conditions.

template <typename Closure, typename ValueType>
void generate_rhs(Closure&& f, ValueType u0, ValueType u1,
{
const auto discretization_points = rhs->get_size()[0];
auto values = rhs->get_values();
Kokkos::View<ValueType*> values_view(values, discretization_points);
Kokkos::parallel_for(
"generate_rhs", discretization_points, KOKKOS_LAMBDA(int i) {
const ValueType h = 1.0 / (discretization_points + 1);
const ValueType xi = ValueType(i + 1) * h;
values_view[i] = -f(xi) * h * h;
if (i == 0) {
values_view[i] += u0;
}
if (i == discretization_points - 1) {
values_view[i] += u1;
}
});
}

Computes the 1-norm of the error given the computed u and the correct solution function correct_u.

template <typename Closure, typename ValueType>
double calculate_error(int discretization_points,
Closure&& correct_u)
{
Kokkos::View<const ValueType*> v_u(u->get_const_values(),
discretization_points);
auto error = 0.0;
Kokkos::parallel_reduce(
"calculate_error", discretization_points,
KOKKOS_LAMBDA(int i, double& lsum) {
const auto h = 1.0 / (discretization_points + 1);
const auto xi = (i + 1) * h;
lsum += Kokkos::Experimental::abs(
(v_u(i) - correct_u(xi)) /
Kokkos::Experimental::abs(correct_u(xi)));
},
error);
return error;
}
int main(int argc, char* argv[])
{
Kokkos::ScopeGuard kokkos(argc, argv);

Some shortcuts

using ValueType = double;
using RealValueType = gko::remove_complex<ValueType>;
using IndexType = int;

Print help message. For details on the kokkos-options see https://kokkos.github.io/kokkos-core-wiki/ProgrammingGuide/Initialization.html#initialization-by-command-line-arguments

if (argc == 2 && (std::string(argv[1]) == "--help")) {
std::cerr << "Usage: " << argv[0]
<< " [discretization_points] [kokkos-options]" << std::endl;
std::exit(-1);
}
const unsigned int discretization_points =
argc >= 2 ? std::atoi(argv[1]) : 100u;

chooses the executor that corresponds to the Kokkos DefaultExecutionSpace

auto exec = []() -> std::shared_ptr<gko::Executor> {
#ifdef KOKKOS_ENABLE_SERIAL
if (std::is_same<Kokkos::DefaultExecutionSpace,
Kokkos::Serial>::value) {
return gko::ReferenceExecutor::create();
}
#endif
#ifdef KOKKOS_ENABLE_OPENMP
if (std::is_same<Kokkos::DefaultExecutionSpace,
Kokkos::OpenMP>::value) {
}
#endif
#ifdef KOKKOS_ENABLE_CUDA
if (std::is_same<Kokkos::DefaultExecutionSpace, Kokkos::Cuda>::value) {
gko::ReferenceExecutor::create());
}
#endif
#ifdef KOKKOS_ENABLE_HIP
if (std::is_same<Kokkos::DefaultExecutionSpace, Kokkos::HIP>::value) {
gko::ReferenceExecutor::create());
}
#endif
}();

problem:

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

initialize vectors

auto rhs = vec::create(exec, gko::dim<2>(discretization_points, 1));
generate_rhs(f, u0, u1, rhs.get());
auto u = vec::create(exec, gko::dim<2>(discretization_points, 1));
for (int i = 0; i < u->get_size()[0]; ++i) {
u->get_values()[i] = 0.0;
}

initialize the stencil matrix

auto A = share(mtx::create(
exec, gko::dim<2>{discretization_points, discretization_points}));
generate_stencil_matrix(A.get());
const RealValueType reduction_factor{1e-7};

Generate solver and solve the system

cg::build()
.with_criteria(
gko::stop::Iteration::build().with_max_iters(discretization_points),
reduction_factor))
.with_preconditioner(bj::build())
.on(exec)
->generate(A)
->apply(rhs, u);
std::cout << "\nSolve complete."
<< "\nThe average relative error is "
<< calculate_error(discretization_points, u.get(), correct_u) /
discretization_points
<< std::endl;
}

Results

Example output:

> ./kokkos-assembly
Solve complete.
The average relative error is 1.05488e-11

The actual error depends on the used hardware.

The plain program

/*******************************<GINKGO LICENSE>******************************
Copyright (c) 2017-2023, 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 <iostream>
#include <map>
#include <string>
#include <omp.h>
#include <Kokkos_Core.hpp>
#include <ginkgo/ginkgo.hpp>
template <typename ValueType, typename IndexType>
void generate_stencil_matrix(gko::matrix::Csr<ValueType, IndexType>* matrix)
{
auto exec = matrix->get_executor();
const auto discretization_points = matrix->get_size()[0];
discretization_points * 3);
Kokkos::View<IndexType*> v_row_idxs(md.get_row_idxs(), md.get_num_elems());
Kokkos::View<IndexType*> v_col_idxs(md.get_col_idxs(), md.get_num_elems());
Kokkos::View<ValueType*> v_values(md.get_values(), md.get_num_elems());
Kokkos::parallel_for(
"generate_stencil_matrix", md.get_num_elems(), KOKKOS_LAMBDA(int i) {
const ValueType coefs[] = {-1, 2, -1};
auto ofs = static_cast<IndexType>((i % 3) - 1);
auto row = static_cast<IndexType>(i / 3);
auto col = row + ofs;
auto mask =
static_cast<IndexType>(0 <= col && col < discretization_points);
v_row_idxs[i] = mask * row;
v_col_idxs[i] = mask * col;
v_values[i] = mask * coefs[ofs + 1];
});
matrix->read(std::move(md));
}
template <typename Closure, typename ValueType>
void generate_rhs(Closure&& f, ValueType u0, ValueType u1,
{
const auto discretization_points = rhs->get_size()[0];
auto values = rhs->get_values();
Kokkos::View<ValueType*> values_view(values, discretization_points);
Kokkos::parallel_for(
"generate_rhs", discretization_points, KOKKOS_LAMBDA(int i) {
const ValueType h = 1.0 / (discretization_points + 1);
const ValueType xi = ValueType(i + 1) * h;
values_view[i] = -f(xi) * h * h;
if (i == 0) {
values_view[i] += u0;
}
if (i == discretization_points - 1) {
values_view[i] += u1;
}
});
}
template <typename Closure, typename ValueType>
double calculate_error(int discretization_points,
Closure&& correct_u)
{
Kokkos::View<const ValueType*> v_u(u->get_const_values(),
discretization_points);
auto error = 0.0;
Kokkos::parallel_reduce(
"calculate_error", discretization_points,
KOKKOS_LAMBDA(int i, double& lsum) {
const auto h = 1.0 / (discretization_points + 1);
const auto xi = (i + 1) * h;
lsum += Kokkos::Experimental::abs(
(v_u(i) - correct_u(xi)) /
Kokkos::Experimental::abs(correct_u(xi)));
},
error);
return error;
}
int main(int argc, char* argv[])
{
Kokkos::ScopeGuard kokkos(argc, argv);
using ValueType = double;
using RealValueType = gko::remove_complex<ValueType>;
using IndexType = int;
if (argc == 2 && (std::string(argv[1]) == "--help")) {
std::cerr << "Usage: " << argv[0]
<< " [discretization_points] [kokkos-options]" << std::endl;
std::exit(-1);
}
const unsigned int discretization_points =
argc >= 2 ? std::atoi(argv[1]) : 100u;
auto exec = []() -> std::shared_ptr<gko::Executor> {
#ifdef KOKKOS_ENABLE_SERIAL
if (std::is_same<Kokkos::DefaultExecutionSpace,
Kokkos::Serial>::value) {
return gko::ReferenceExecutor::create();
}
#endif
#ifdef KOKKOS_ENABLE_OPENMP
if (std::is_same<Kokkos::DefaultExecutionSpace,
Kokkos::OpenMP>::value) {
}
#endif
#ifdef KOKKOS_ENABLE_CUDA
if (std::is_same<Kokkos::DefaultExecutionSpace, Kokkos::Cuda>::value) {
gko::ReferenceExecutor::create());
}
#endif
#ifdef KOKKOS_ENABLE_HIP
if (std::is_same<Kokkos::DefaultExecutionSpace, Kokkos::HIP>::value) {
gko::ReferenceExecutor::create());
}
#endif
}();
auto correct_u = [] KOKKOS_FUNCTION(ValueType x) { return x * x * x; };
auto f = [] KOKKOS_FUNCTION(ValueType x) { return ValueType{6} * x; };
auto u0 = correct_u(0);
auto u1 = correct_u(1);
auto rhs = vec::create(exec, gko::dim<2>(discretization_points, 1));
generate_rhs(f, u0, u1, rhs.get());
auto u = vec::create(exec, gko::dim<2>(discretization_points, 1));
for (int i = 0; i < u->get_size()[0]; ++i) {
u->get_values()[i] = 0.0;
}
auto A = share(mtx::create(
exec, gko::dim<2>{discretization_points, discretization_points}));
generate_stencil_matrix(A.get());
const RealValueType reduction_factor{1e-7};
cg::build()
.with_criteria(
gko::stop::Iteration::build().with_max_iters(discretization_points),
reduction_factor))
.with_preconditioner(bj::build())
.on(exec)
->generate(A)
->apply(rhs, u);
std::cout << "\nSolve complete."
<< "\nThe average relative error is "
<< calculate_error(discretization_points, u.get(), correct_u) /
discretization_points
<< 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:54
gko::matrix::Dense
Dense is a matrix format which explicitly stores all values of the matrix.
Definition: dense_cache.hpp:48
gko::device_matrix_data::sum_duplicates
void sum_duplicates()
Sums up all duplicate entries pointing to the same non-zero location.
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::preconditioner::Jacobi
A block-Jacobi preconditioner is a block-diagonal linear operator, obtained by inverting the diagonal...
Definition: jacobi.hpp:213
gko::stop::ResidualNorm
The ResidualNorm class is a stopping criterion which stops the iteration process when the actual resi...
Definition: residual_norm.hpp:138
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:74
gko::share
detail::shared_type< OwningPointer > share(OwningPointer &&p)
Marks the object pointed to by p as shared.
Definition: utils_helper.hpp:254
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::solver::initial_guess_mode::rhs
the input is right hand side
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:1373
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:865
gko::PolymorphicObject::get_executor
std::shared_ptr< const Executor > get_executor() const noexcept
Returns the Executor of the object.
Definition: polymorphic_object.hpp:263
gko::LinOp::get_size
const dim< 2 > & get_size() const noexcept
Returns the size of the operator.
Definition: lin_op.hpp:239
gko::remove_complex
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:354
gko::device_matrix_data
This type is a device-side equivalent to matrix_data.
Definition: device_matrix_data.hpp:63
gko::matrix::Csr::read
void read(const mat_data &data) override
Reads a matrix from a matrix_data structure.