Ginkgo  Generated from pipelines/224724463 branch based on develop. Ginkgo version 1.3.0
A numerical linear algebra library targeting many-core architectures
The minimal-cuda-solver program

The minimal CUDA 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 is a minimal example that solves a system with Ginkgo. The matrix, right hand side and initial guess are read from standard input, and the result is written to standard output. The system matrix is stored in CSR format, and the system solved using the CG method, preconditioned with the block-Jacobi preconditioner. All computations are done on the GPU.

The easiest way to use the example data from the data/ folder is to concatenate the matrix, the right hand side and the initial solution (in that exact order), and pipe the result to the minimal_solver_cuda executable:

cat data/A.mtx data/b.mtx data/x0.mtx | ./minimal-cuda-solver

About the example

The commented program

#include <ginkgo/ginkgo.hpp>
#include <iostream>
int main()
{

Instantiate a CUDA executor

Read data

auto A = gko::read<gko::matrix::Csr<>>(std::cin, gpu);
auto b = gko::read<gko::matrix::Dense<>>(std::cin, gpu);
auto x = gko::read<gko::matrix::Dense<>>(std::cin, gpu);

Create the solver

auto solver =
.with_preconditioner(gko::preconditioner::Jacobi<>::build().on(gpu))
.with_criteria(
gko::stop::Iteration::build().with_max_iters(20u).on(gpu),
.with_reduction_factor(1e-15)
.on(gpu))
.on(gpu);

Solve system

solver->generate(give(A))->apply(lend(b), lend(x));

Write result

write(std::cout, lend(x));
}

Results

The following is the expected result when using the data contained in the folder data as input:

%%MatrixMarket matrix array real general
19 1
0.252218
0.108645
0.0662811
0.0630433
0.0384088
0.0396536
0.0402648
0.0338935
0.0193098
0.0234653
0.0211499
0.0196413
0.0199151
0.0181674
0.0162722
0.0150714
0.0107016
0.0121141
0.0123025

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>
int main()
{
auto A = gko::read<gko::matrix::Csr<>>(std::cin, gpu);
auto b = gko::read<gko::matrix::Dense<>>(std::cin, gpu);
auto x = gko::read<gko::matrix::Dense<>>(std::cin, gpu);
auto solver =
.with_preconditioner(gko::preconditioner::Jacobi<>::build().on(gpu))
.with_criteria(
gko::stop::Iteration::build().with_max_iters(20u).on(gpu),
.with_reduction_factor(1e-15)
.on(gpu))
.on(gpu);
solver->generate(give(A))->apply(lend(b), lend(x));
write(std::cout, lend(x));
}
gko::real
constexpr std::enable_if_t<!is_complex_s< T >::value, T > real(const T &x)
Returns the real part of the object.
Definition: math.hpp:817
gko::layout_type::array
The matrix should be written as dense matrix in column-major order.
gko::give
std::remove_reference< OwningPointer >::type && give(OwningPointer &&p)
Marks that the object pointed to by p can be given to the callee.
Definition: utils.hpp:231
gko::CudaExecutor::create
static std::shared_ptr< CudaExecutor > create(int device_id, std::shared_ptr< Executor > master, bool device_reset=false)
Creates a new CudaExecutor.
gko::OmpExecutor::create
static std::shared_ptr< OmpExecutor > create()
Creates a new OmpExecutor.
Definition: executor.hpp:909
gko::stop::ResidualNormReduction
The ResidualNormReduction class is a stopping criterion which stops the iteration process when the re...
Definition: residual_norm.hpp:113
gko::write
void write(StreamType &&os, MatrixType *matrix, layout_type layout=layout_type::array)
Reads a matrix stored in matrix market format from an input stream.
Definition: mtx_io.hpp:134
gko::preconditioner::Jacobi
A block-Jacobi preconditioner is a block-diagonal linear operator, obtained by inverting the diagonal...
Definition: jacobi.hpp:207
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::solver::Cg
CG or the conjugate gradient method is an iterative type Krylov subspace method which is suitable for...
Definition: cg.hpp:72