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

The preconditioned solver example..

This example depends on simple-solver.

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

Introduction

About the example

The commented program

#include <ginkgo/ginkgo.hpp>
#include <fstream>
#include <iostream>
#include <string>
int main(int argc, char *argv[])
{

Some shortcuts

Print version information

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

Figure out where to run the code

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);
}

Read data

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);

Create solver factory

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

Add preconditioner, these 2 lines are the only difference from the simple solver example

.with_preconditioner(bj::build().with_max_block_size(8u).on(exec))
.on(exec);

Create solver

auto solver = solver_gen->generate(A);

Solve system

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

Print solution

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

Calculate residual

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

This is the expected output:

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
9.08137e-16

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-20)
.on(exec))
.with_preconditioner(bj::build().with_max_block_size(8u).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::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
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