The poisson solver example.
This example depends on simple-solver.
Introduction
This example solves a 1D Poisson equation:
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:
For an equidistant grid with K "inner" discretization points and step size , the formula produces a system of linear equations
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 is set to (making the solution ), but that can be changed in the main
function.
The intention of the example is to show how Ginkgo can be used to build an application solving a real-world problem, which includes a solution of a large, sparse linear system as a component.
About the example
The commented program
row_ptrs[i + 1] = pos;
}
}
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();
const ValueType h = 1.0 / static_cast< ValueType> (discretization_points + 1);
const auto xi = static_cast< ValueType> (i + 1) * h;
values[i] = -f(xi) * h * h;
}
values[0] += u0;
values[discretization_points - 1] += u1;
}
Dense is a matrix format which explicitly stores all values of the matrix.
Definition dense.hpp:117
std::size_t size_type
Integral type used for allocation quantities.
Definition types.hpp:89
Prints the solution u
.
template <typename Closure, typename ValueType>
void print_solution(ValueType u0, ValueType u1,
{
std::cout << u0 << '\n' ;
for (
int i = 0; i < u->
get_size ()[0]; ++i) {
}
std::cout << u1 << std::endl;
}
const dim< 2 > & get_size() const noexcept
Returns the size of the operator.
Definition lin_op.hpp:210
const value_type * get_const_values() const noexcept
Returns a pointer to the array of values of the matrix.
Definition dense.hpp:860
Computes the 1-norm of the error given the computed u
and the correct solution function correct_u
.
template <typename Closure, typename ValueType>
Closure correct_u)
{
const ValueType h = 1.0 / static_cast< ValueType> (discretization_points + 1);
for (int i = 0; i < discretization_points; ++i) {
using std::abs;
const auto xi = static_cast< ValueType> (i + 1) * h;
error +=
}
return error;
}
int main(int argc, char * argv[])
{
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 ValueType = double;
using IndexType = int;
CSR is a matrix format which stores only the nonzero coefficients by compressing each row of the matr...
Definition csr.hpp:121
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
Print version information
if (argc == 2 && (std::string(argv[1]) == "--help" )) {
std::cerr << "Usage: " << argv[0]
<< " [executor] [DISCRETIZATION_POINTS]" << std::endl;
std::exit(-1);
}
static const version_info & get()
Returns an instance of version_info.
Definition version.hpp:139
Get number of discretization points
const unsigned int discretization_points =
argc >= 3 ? std::atoi(argv[2]) : 100;
Get the executor string
const auto executor_string = argc >= 2 ? argv[1] : "reference" ;
Figure out where to run the code
std::map<std::string, std::function<std::shared_ptr<gko::Executor>()>>
exec_map{
{"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)();
executor used by the application
const auto app_exec = exec->get_master();
Set up the problem: define the exact solution, the right hand side and the Dirichlet boundary condition.
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);
initialize matrix and vectors
auto matrix = mtx::create(app_exec,
gko::dim<2> (discretization_points),
3 * discretization_points - 2);
generate_stencil_matrix(matrix.get());
auto rhs = vec::create(app_exec,
gko::dim<2> (discretization_points, 1));
generate_rhs(f, u0, u1, rhs.get());
auto u = vec::create(app_exec,
gko::dim<2> (discretization_points, 1));
for (
int i = 0; i < u->
get_size ()[0]; ++i) {
}
value_type * get_values() noexcept
Returns a pointer to the array of values of the matrix.
Definition dense.hpp:851
A type representing the dimensions of a multidimensional object.
Definition dim.hpp:26
Generate solver and solve the system
cg::build()
.with_criteria(
gko::stop::Iteration::build().with_max_iters(discretization_points),
gko::stop::ResidualNorm<ValueType>::build().with_reduction_factor(
reduction_factor))
.with_preconditioner(bj::build())
.on(exec)
->generate(clone(exec, matrix))
->apply(rhs, u);
Uncomment to print the solution print_solution<ValueType>(u0, u1, u.get());
std::cout << "Solve complete.\nThe average relative error is "
<< calculate_error(discretization_points, u.get(), correct_u) /
discretization_points)
<< std::endl;
}
Results
This is the expected output:
Solve complete.
The average relative error is 2.52236e-11
Comments about programming and debugging
The plain program
#include <iostream>
#include <map>
#include <string>
#include <vector>
#include <ginkgo/ginkgo.hpp>
template <typename ValueType, typename IndexType>
{
const auto discretization_points = matrix->
get_size ()[0];
int pos = 0;
const ValueType coefs[] = {-1, 2, -1};
row_ptrs[0] = pos;
for (int 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>
void generate_rhs(Closure f, ValueType u0, ValueType u1,
{
const auto discretization_points =
rhs ->get_size()[0];
auto values =
rhs ->get_values();
const ValueType h = 1.0 / static_cast< ValueType> (discretization_points + 1);
const auto xi = static_cast< ValueType> (i + 1) * h;
values[i] = -f(xi) * h * h;
}
values[0] += u0;
values[discretization_points - 1] += u1;
}
template <typename Closure, typename ValueType>
void print_solution(ValueType u0, ValueType u1,
{
std::cout << u0 << '\n' ;
for (
int i = 0; i < u->
get_size ()[0]; ++i) {
}
std::cout << u1 << std::endl;
}
template <typename Closure, typename ValueType>
Closure correct_u)
{
const ValueType h = 1.0 / static_cast< ValueType> (discretization_points + 1);
for (int i = 0; i < discretization_points; ++i) {
using std::abs;
const auto xi = static_cast< ValueType> (i + 1) * h;
error +=
}
return error;
}
int main(int argc, char * argv[])
{
using ValueType = double;
using IndexType = int;
if (argc == 2 && (std::string(argv[1]) == "--help" )) {
std::cerr << "Usage: " << argv[0]
<< " [executor] [DISCRETIZATION_POINTS]" << std::endl;
std::exit(-1);
}
const unsigned int discretization_points =
argc >= 3 ? std::atoi(argv[2]) : 100;
const auto executor_string = argc >= 2 ? argv[1] : "reference" ;
std::map<std::string, std::function<std::shared_ptr<gko::Executor>()>>
exec_map{
{"cuda" ,
[] {
}},
{"hip" ,
[] {
}},
{"dpcpp" ,
[] {
}},
{"reference" , [] { return gko::ReferenceExecutor::create(); }}};
const auto exec = exec_map.at(executor_string)();
const auto app_exec = exec->get_master();
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);
auto matrix = mtx::create(app_exec,
gko::dim<2> (discretization_points),
3 * discretization_points - 2);
generate_stencil_matrix(matrix.get());
auto rhs = vec::create(app_exec,
gko::dim<2> (discretization_points, 1));
generate_rhs(f, u0, u1,
rhs .get());
auto u = vec::create(app_exec,
gko::dim<2> (discretization_points, 1));
for (
int i = 0; i < u->
get_size ()[0]; ++i) {
}
cg::build()
.with_criteria(
gko::stop::Iteration::build().with_max_iters(discretization_points),
gko::stop::ResidualNorm<ValueType>::build().with_reduction_factor(
reduction_factor))
.with_preconditioner(bj::build())
.on(exec)
->generate(
clone (exec, matrix))
->apply(rhs, u);
std::cout << "Solve complete.\nThe average relative error is "
<< calculate_error(discretization_points, u.get(), correct_u) /
discretization_points)
<< std::endl;
}
index_type * get_row_ptrs() noexcept
Returns the row pointers of the matrix.
Definition csr.hpp:905
value_type * get_values() noexcept
Returns the values of the matrix.
Definition csr.hpp:867
index_type * get_col_idxs() noexcept
Returns the column indexes of the matrix.
Definition csr.hpp:886
@ 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
detail::cloned_type< Pointer > clone(const Pointer &p)
Creates a unique clone of the object pointed to by p.
Definition utils_helper.hpp:173