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

The simple solver example..

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

Introduction

This simple solver example should help you get started with Ginkgo. This example is meant for you to understand how Ginkgo works and how you can solve a simple linear system with Ginkgo. We encourage you to play with the code, change the parameters and see what is best suited for your purposes.

About the example

Each example has the following sections:

  1. Introduction:This gives an overview of the example and mentions any interesting aspects in the example that might help the reader.
  2. The commented program: This section is intended for you to understand the details of the example so that you can play with it and understand Ginkgo and its features better.
  3. Results: This section shows the results of the code when run. Though the results may not be completely the same, you can expect the behaviour to be similar.
  4. The plain program: This is the complete code without any comments to have an complete overview of the code.

The commented program

Include files

This is the main ginkgo header file.

#include <ginkgo/ginkgo.hpp>

Add the fstream header to read from data from files.

#include <fstream>

Add the C++ iostream header to output information to the console.

#include <iostream>

Add the string manipulation header to handle strings.

#include <string>
int main(int argc, char *argv[])
{

Use some shortcuts. In Ginkgo, vectors are seen as a gko::matrix::Dense with one column/one row. The advantage of this concept is that using multiple vectors is a now a natural extension of adding columns/rows are necessary.

The gko::matrix::Csr class is used here, but any other matrix class such as gko::matrix::Coo, gko::matrix::Hybrid, gko::matrix::Ell or gko::matrix::Sellp could also be used.

using mtx = gko::matrix::Csr<>;

The gko::solver::Cg is used here, but any other solver class can also be used.

using cg = gko::solver::Cg<>;

Print the ginkgo version information.

std::cout << gko::version_info::get() << std::endl;

Where do you want to run your solver ?

The gko::Executor class is one of the cornerstones of Ginkgo. Currently, we have support for an gko::OmpExecutor, which uses OpenMP multi-threading in most of its kernels, a gko::ReferenceExecutor, a single threaded specialization of the OpenMP executor and a gko::CudaExecutor which runs the code on a NVIDIA GPU if available.

Note
With the help of C++, you see that you only ever need to change the executor and all the other functions/ routines within Ginkgo should automatically work and run on the executor with any other changes.
std::shared_ptr<gko::Executor> exec;
if (argc == 1 || std::string(argv[1]) == "reference") {
exec = gko::ReferenceExecutor::create();
} else if (argc == 2 && std::string(argv[1]) == "omp") {
} else if (argc == 2 && std::string(argv[1]) == "cuda" &&
} else if (argc == 2 && std::string(argv[1]) == "hip" &&
} else {
std::cerr << "Usage: " << argv[0] << " [executor]" << std::endl;
std::exit(-1);
}

Reading your data and transfer to the proper device.

Read the matrix, right hand side and the initial solution using the read function.

Note
Ginkgo uses C++ smart pointers to automatically manage memory. To this end, we use our own object ownership transfer functions that under the hood call the required smart pointer functions to manage object ownership. The gko::share , gko::give and gko::lend are the functions that you would need to use.
auto A = share(gko::read<mtx>(std::ifstream("data/A.mtx"), exec));
auto b = gko::read<vec>(std::ifstream("data/b.mtx"), exec);
auto x = gko::read<vec>(std::ifstream("data/x0.mtx"), exec);

Creating the solver

Generate the gko::solver factory. Ginkgo uses the concept of Factories to build solvers with certain properties. Observe the Fluent interface used here. Here a cg solver is generated with a stopping criteria of maximum iterations of 20 and a residual norm reduction of 1e-15. You also observe that the stopping criteria(gko::stop) are also generated from factories using their build methods. You need to specify the executors which each of the object needs to be built on.

auto solver_gen =
cg::build()
.with_criteria(
gko::stop::Iteration::build().with_max_iters(20u).on(exec),
.with_reduction_factor(1e-15)
.on(exec))
.on(exec);

Generate the solver from the matrix. The solver factory built in the previous step takes a "matrix"(a gko::LinOp to be more general) as an input. In this case we provide it with a full matrix that we previously read, but as the solver only effectively uses the apply() method within the provided "matrix" object, you can effectively create a gko::LinOp class with your own apply implementation to accomplish more tasks. We will see an example of how this can be done in the custom-matrix-format example

auto solver = solver_gen->generate(A);

Finally, solve the system. The solver, being a gko::LinOp, can be applied to a right hand side, b to obtain the solution, x.

solver->apply(lend(b), lend(x));

Print the solution to the command line.

std::cout << "Solution (x): \n";
write(std::cout, lend(x));

To measure if your solution has actually converged, you can measure the error of the solution. one, neg_one are objects that represent the numbers which allow for a uniform interface when computing on any device. To compute the residual, all you need to do is call the apply method, which in this case is an spmv and equivalent to the LAPACK z_spmv routine. Finally, you compute the euclidean 2-norm with the compute_norm2 function.

auto one = gko::initialize<vec>({1.0}, exec);
auto neg_one = gko::initialize<vec>({-1.0}, exec);
auto res = gko::initialize<vec>({0.0}, exec);
A->apply(lend(one), lend(x), lend(neg_one), lend(b));
b->compute_norm2(lend(res));
std::cout << "Residual norm sqrt(r^T r): \n";
write(std::cout, lend(res));
}

Results

The following is the expected result:

Solution (x):
%%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
Residual norm sqrt(r^T r):
%%MatrixMarket matrix array real general
1 1
2.10788e-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 <fstream>
#include <iostream>
#include <string>
int main(int argc, char *argv[])
{
using vec = gko::matrix::Dense<>;
using mtx = gko::matrix::Csr<>;
using cg = gko::solver::Cg<>;
std::cout << gko::version_info::get() << std::endl;
std::shared_ptr<gko::Executor> exec;
if (argc == 1 || std::string(argv[1]) == "reference") {
exec = gko::ReferenceExecutor::create();
} else if (argc == 2 && std::string(argv[1]) == "omp") {
} else if (argc == 2 && std::string(argv[1]) == "cuda" &&
} else if (argc == 2 && std::string(argv[1]) == "hip" &&
} else {
std::cerr << "Usage: " << argv[0] << " [executor]" << std::endl;
std::exit(-1);
}
auto A = share(gko::read<mtx>(std::ifstream("data/A.mtx"), exec));
auto b = gko::read<vec>(std::ifstream("data/b.mtx"), exec);
auto x = gko::read<vec>(std::ifstream("data/x0.mtx"), exec);
auto solver_gen =
cg::build()
.with_criteria(
gko::stop::Iteration::build().with_max_iters(20u).on(exec),
.with_reduction_factor(1e-15)
.on(exec))
.on(exec);
auto solver = solver_gen->generate(A);
solver->apply(lend(b), lend(x));
std::cout << "Solution (x): \n";
write(std::cout, lend(x));
auto one = gko::initialize<vec>({1.0}, exec);
auto neg_one = gko::initialize<vec>({-1.0}, exec);
auto res = gko::initialize<vec>({0.0}, exec);
A->apply(lend(one), lend(x), lend(neg_one), lend(b));
b->compute_norm2(lend(res));
std::cout << "Residual norm sqrt(r^T r): \n";
write(std::cout, lend(res));
}
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::CudaExecutor::get_num_devices
static int get_num_devices()
Get the number of devices present on the system.
gko::layout_type::array
The matrix should be written as dense matrix in column-major order.
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::HipExecutor::get_num_devices
static int get_num_devices()
Get the number of devices present on the system.
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::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::version_info::get
static const version_info & get()
Returns an instance of version_info.
Definition: version.hpp:168
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
gko::share
detail::shared_type< OwningPointer > share(OwningPointer &&p)
Marks the object pointed to by p as shared.
Definition: utils.hpp:210
gko::real
constexpr xstd::enable_if_t<!is_complex_s< T >::value, T > real(const T &x)
Returns the real part of the object.
Definition: math.hpp:630
gko::one
constexpr T one()
Returns the multiplicative identity for T.
Definition: math.hpp:534