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

The simple solver with a custom logger example..

This example depends on simple-solver, simple-solver-logging, minimal-cuda-solver.

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

Introduction

The custom-logger example shows how Ginkgo's API can be leveraged to implement application-specific callbacks for Ginkgo's events. This is the most basic way of extending Ginkgo and a good first step for any application developer who wants to adapt Ginkgo to his specific needs.

Ginkgo's gko::log::Logger abstraction provides hooks to the events that happen during the library execution. These hooks concern any low-level event such as memory allocations, deallocations, copies and kernel launches up to high-level events such as linear operator applications and completion of solver iterations.

In this example, a simple logger is implemented to track the solver's recurrent residual norm and compute the true residual norm. At the end of the solver execution, a comparison table is shown on-screen.

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 map header for storing the executor map.

#include <map>

Add the C++ iomanip header to prettify the output.

#include <iomanip>

Add formatting flag modification capabilities.

#include <ios>

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

#include <iostream>

Add the string manipulation header to handle strings.

#include <string>

Add the vector header for storing the logger's data

#include <vector>

Utility function which returns the first element (position [0, 0]) from a given gko::matrix::Dense matrix / vector.

template <typename ValueType>
ValueType get_first_element(const gko::matrix::Dense<ValueType>* mtx)
{

Copy the matrix / vector to the host device before accessing the value in case it is stored in a GPU.

return mtx->get_executor()->copy_val_to_host(mtx->get_const_values());
}

Utility function which computes the norm of a Ginkgo gko::matrix::Dense vector.

template <typename ValueType>
{

Get the executor of the vector

auto exec = b->get_executor();

Initialize a result scalar containing the value 0.0.

auto b_norm =
gko::initialize<gko::matrix::Dense<gko::remove_complex<ValueType>>>(
{0.0}, exec);

Use the dense compute_norm2 function to compute the norm.

b->compute_norm2(b_norm);

Use the other utility function to return the norm contained in b_norm

return get_first_element(b_norm.get());
}

Custom logger class which intercepts the residual norm scalar and solution vector in order to print a table of real vs recurrent (internal to the solvers) residual norms.

template <typename ValueType>
struct ResidualLogger : gko::log::Logger {
using RealValueType = gko::remove_complex<ValueType>;

Output the logger's data in a table format

void write() const
{

Print a header for the table

std::cout << "Recurrent vs true vs implicit residual norm:"
<< std::endl;
std::cout << '|' << std::setw(10) << "Iteration" << '|' << std::setw(25)
<< "Recurrent Residual Norm" << '|' << std::setw(25)
<< "True Residual Norm" << '|' << std::setw(25)
<< "Implicit Residual Norm" << '|' << std::endl;

Print a separation line. Note that for creating 10 characters std::setw() should be set to 11.

std::cout << '|' << std::setfill('-') << std::setw(11) << '|'
<< std::setw(26) << '|' << std::setw(26) << '|'
<< std::setw(26) << '|' << std::setfill(' ') << std::endl;

Print the data one by one in the form

std::cout << std::scientific;
for (std::size_t i = 0; i < iterations.size(); i++) {
std::cout << '|' << std::setw(10) << iterations[i] << '|'
<< std::setw(25) << recurrent_norms[i] << '|'
<< std::setw(25) << real_norms[i] << '|' << std::setw(25)
<< implicit_norms[i] << '|' << std::endl;
}

std::defaultfloat could be used here but some compilers do not support it properly, e.g. the Intel compiler

std::cout.unsetf(std::ios_base::floatfield);

Print a separation line

std::cout << '|' << std::setfill('-') << std::setw(11) << '|'
<< std::setw(26) << '|' << std::setw(26) << '|'
<< std::setw(26) << '|' << std::setfill(' ') << std::endl;
}
using gko_real_dense = gko::matrix::Dense<RealValueType>;

Customize the logging hook which is called everytime an iteration is completed

void on_iteration_complete(const gko::LinOp* solver, const gko::LinOp* b,
const gko::LinOp* solution,
const gko::size_type& iteration,
const gko::LinOp* residual,
const gko::LinOp* residual_norm,
const gko::LinOp* implicit_sq_residual_norm,
bool) const override
{

If the solver shares a residual norm, log its value

if (residual_norm) {
auto dense_norm = gko::as<gko_real_dense>(residual_norm);

Add the norm to the recurrent_norms vector

recurrent_norms.push_back(get_first_element(dense_norm));

Otherwise, use the recurrent residual vector

} else {
auto dense_residual = gko::as<gko_dense>(residual);

Compute the residual vector's norm

auto norm = compute_norm(dense_residual);

Add the computed norm to the recurrent_norms vector

recurrent_norms.push_back(norm);
}

If the solver shares the current solution vector

if (solution) {

Extract the matrix from the solver

auto matrix = gko::as<gko::solver::detail::SolverBaseLinOp>(solver)
->get_system_matrix();

Store the matrix's executor

auto exec = matrix->get_executor();

Create a scalar containing the value 1.0

auto one = gko::initialize<gko_dense>({1.0}, exec);

Create a scalar containing the value -1.0

auto neg_one = gko::initialize<gko_dense>({-1.0}, exec);

Instantiate a temporary result variable

auto res = gko::as<gko_dense>(gko::clone(b));

Compute the real residual vector by calling apply on the system matrix

matrix->apply(one, solution, neg_one, res);

Compute the norm of the residual vector and add it to the real_norms vector

real_norms.push_back(compute_norm(res.get()));
} else {

Add to the real_norms vector the value -1.0 if it could not be computed

real_norms.push_back(-1.0);
}
if (implicit_sq_residual_norm) {
auto dense_norm =
gko::as<gko_real_dense>(implicit_sq_residual_norm);

Add the norm to the implicit_norms vector

implicit_norms.push_back(std::sqrt(get_first_element(dense_norm)));
} else {

Add to the implicit_norms vector the value -1.0 if it could not be computed

implicit_norms.push_back(-1.0);
}

Add the current iteration number to the iterations vector

iterations.push_back(iteration);
}

Construct the logger

ResidualLogger()
: gko::log::Logger(gko::log::Logger::iteration_complete_mask)
{}
private:

Vector which stores all the recurrent residual norms

mutable std::vector<RealValueType> recurrent_norms{};

Vector which stores all the real residual norms

mutable std::vector<RealValueType> real_norms{};

Vector which stores all the implicit residual norms

mutable std::vector<RealValueType> implicit_norms{};

Vector which stores all the iteration numbers

mutable std::vector<std::size_t> iterations{};
};
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.

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

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.

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

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.
if (argc == 2 && (std::string(argv[1]) == "--help")) {
std::cerr << "Usage: " << argv[0] << " [executor]" << std::endl;
std::exit(-1);
}
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{
{"omp", [] { return gko::OmpExecutor::create(); }},
{"cuda",
[] {
}},
{"hip",
[] {
}},
{"dpcpp",
[] {
}},
{"reference", [] { return gko::ReferenceExecutor::create(); }}};

executor where Ginkgo will perform the computation

const auto exec = exec_map.at(executor_string)(); // throws if not valid

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. gko::share and gko::give 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);
const RealValueType reduction_factor = 1e-7;

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),
.with_reduction_factor(reduction_factor))
.on(exec);

Instantiate a ResidualLogger logger.

auto logger = std::make_shared<ResidualLogger<ValueType>>();

Add the previously created logger to the solver factory. The logger will be automatically propagated to all solvers created from this factory.

solver_gen->add_logger(logger);

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(b, x);

Print the solution to the command line.

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

Print the table of the residuals obtained from the logger

logger->write();

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<real_vec>({0.0}, exec);
A->apply(one, x, neg_one, b);
b->compute_norm2(res);
std::cout << "Residual norm sqrt(r^T r):\n";
write(std::cout, 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
Recurrent vs true vs implicit residual norm:
| Iteration| Recurrent Residual Norm| True Residual Norm| Implicit Residual Norm|
|----------|-------------------------|-------------------------|-------------------------|
| 0| 4.358899e+00| 4.358899e+00| 4.358899e+00|
| 1| 2.304548e+00| 2.304548e+00| 2.304548e+00|
| 2| 1.467706e+00| 1.467706e+00| 1.467706e+00|
| 3| 9.848751e-01| 9.848751e-01| 9.848751e-01|
| 4| 7.418330e-01| 7.418330e-01| 7.418330e-01|
| 5| 5.136231e-01| 5.136231e-01| 5.136231e-01|
| 6| 3.841650e-01| 3.841650e-01| 3.841650e-01|
| 7| 3.164394e-01| 3.164394e-01| 3.164394e-01|
| 8| 2.277088e-01| 2.277088e-01| 2.277088e-01|
| 9| 1.703121e-01| 1.703121e-01| 1.703121e-01|
| 10| 9.737220e-02| 9.737220e-02| 9.737220e-02|
| 11| 6.168306e-02| 6.168306e-02| 6.168306e-02|
| 12| 4.541231e-02| 4.541231e-02| 4.541231e-02|
| 13| 3.195304e-02| 3.195304e-02| 3.195304e-02|
| 14| 1.616058e-02| 1.616058e-02| 1.616058e-02|
| 15| 6.570152e-03| 6.570152e-03| 6.570152e-03|
| 16| 2.643669e-03| 2.643669e-03| 2.643669e-03|
| 17| 8.588089e-04| 8.588089e-04| 8.588089e-04|
| 18| 2.864613e-04| 2.864613e-04| 2.864613e-04|
| 19| 1.641952e-15| 2.107881e-15| 1.641952e-15|
|----------|-------------------------|-------------------------|-------------------------|
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-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 <ginkgo/ginkgo.hpp>
#include <fstream>
#include <map>
#include <iomanip>
#include <ios>
#include <iostream>
#include <string>
#include <vector>
template <typename ValueType>
ValueType get_first_element(const gko::matrix::Dense<ValueType>* mtx)
{
return mtx->get_executor()->copy_val_to_host(mtx->get_const_values());
}
template <typename ValueType>
{
auto exec = b->get_executor();
auto b_norm =
gko::initialize<gko::matrix::Dense<gko::remove_complex<ValueType>>>(
{0.0}, exec);
b->compute_norm2(b_norm);
return get_first_element(b_norm.get());
}
template <typename ValueType>
struct ResidualLogger : gko::log::Logger {
using RealValueType = gko::remove_complex<ValueType>;
void write() const
{
std::cout << "Recurrent vs true vs implicit residual norm:"
<< std::endl;
std::cout << '|' << std::setw(10) << "Iteration" << '|' << std::setw(25)
<< "Recurrent Residual Norm" << '|' << std::setw(25)
<< "True Residual Norm" << '|' << std::setw(25)
<< "Implicit Residual Norm" << '|' << std::endl;
std::cout << '|' << std::setfill('-') << std::setw(11) << '|'
<< std::setw(26) << '|' << std::setw(26) << '|'
<< std::setw(26) << '|' << std::setfill(' ') << std::endl;
std::cout << std::scientific;
for (std::size_t i = 0; i < iterations.size(); i++) {
std::cout << '|' << std::setw(10) << iterations[i] << '|'
<< std::setw(25) << recurrent_norms[i] << '|'
<< std::setw(25) << real_norms[i] << '|' << std::setw(25)
<< implicit_norms[i] << '|' << std::endl;
}
std::cout.unsetf(std::ios_base::floatfield);
std::cout << '|' << std::setfill('-') << std::setw(11) << '|'
<< std::setw(26) << '|' << std::setw(26) << '|'
<< std::setw(26) << '|' << std::setfill(' ') << std::endl;
}
using gko_dense = gko::matrix::Dense<ValueType>;
using gko_real_dense = gko::matrix::Dense<RealValueType>;
void on_iteration_complete(const gko::LinOp* solver, const gko::LinOp* b,
const gko::LinOp* solution,
const gko::size_type& iteration,
const gko::LinOp* residual,
const gko::LinOp* residual_norm,
const gko::LinOp* implicit_sq_residual_norm,
bool) const override
{
if (residual_norm) {
auto dense_norm = gko::as<gko_real_dense>(residual_norm);
recurrent_norms.push_back(get_first_element(dense_norm));
} else {
auto dense_residual = gko::as<gko_dense>(residual);
auto norm = compute_norm(dense_residual);
recurrent_norms.push_back(norm);
}
if (solution) {
auto matrix = gko::as<gko::solver::detail::SolverBaseLinOp>(solver)
->get_system_matrix();
auto exec = matrix->get_executor();
auto one = gko::initialize<gko_dense>({1.0}, exec);
auto neg_one = gko::initialize<gko_dense>({-1.0}, exec);
auto res = gko::as<gko_dense>(gko::clone(b));
matrix->apply(one, solution, neg_one, res);
real_norms.push_back(compute_norm(res.get()));
} else {
real_norms.push_back(-1.0);
}
if (implicit_sq_residual_norm) {
auto dense_norm =
gko::as<gko_real_dense>(implicit_sq_residual_norm);
implicit_norms.push_back(std::sqrt(get_first_element(dense_norm)));
} else {
implicit_norms.push_back(-1.0);
}
iterations.push_back(iteration);
}
ResidualLogger()
: gko::log::Logger(gko::log::Logger::iteration_complete_mask)
{}
private:
mutable std::vector<RealValueType> recurrent_norms{};
mutable std::vector<RealValueType> real_norms{};
mutable std::vector<RealValueType> implicit_norms{};
mutable std::vector<std::size_t> iterations{};
};
int main(int argc, char* argv[])
{
using ValueType = double;
using RealValueType = gko::remove_complex<ValueType>;
using IndexType = int;
std::cout << gko::version_info::get() << std::endl;
if (argc == 2 && (std::string(argv[1]) == "--help")) {
std::cerr << "Usage: " << argv[0] << " [executor]" << std::endl;
std::exit(-1);
}
const auto executor_string = argc >= 2 ? argv[1] : "reference";
std::map<std::string, std::function<std::shared_ptr<gko::Executor>()>>
exec_map{
{"omp", [] { return gko::OmpExecutor::create(); }},
{"cuda",
[] {
}},
{"hip",
[] {
}},
{"dpcpp",
[] {
}},
{"reference", [] { return gko::ReferenceExecutor::create(); }}};
const auto exec = exec_map.at(executor_string)(); // throws if not valid
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);
const RealValueType reduction_factor = 1e-7;
auto solver_gen =
cg::build()
.with_criteria(gko::stop::Iteration::build().with_max_iters(20u),
.with_reduction_factor(reduction_factor))
.on(exec);
auto logger = std::make_shared<ResidualLogger<ValueType>>();
solver_gen->add_logger(logger);
auto solver = solver_gen->generate(A);
solver->apply(b, x);
std::cout << "Solution (x):\n";
write(std::cout, x);
logger->write();
auto one = gko::initialize<vec>({1.0}, exec);
auto neg_one = gko::initialize<vec>({-1.0}, exec);
auto res = gko::initialize<real_vec>({0.0}, exec);
A->apply(one, x, neg_one, b);
b->compute_norm2(res);
std::cout << "Residual norm sqrt(r^T r):\n";
write(std::cout, res);
}
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::log::profile_event_category::solver
Solver events.
gko::layout_type::array
The matrix should be written as dense matrix in column-major order.
gko::LinOp
Definition: lin_op.hpp:146
gko::matrix::Dense
Dense is a matrix format which explicitly stores all values of the matrix.
Definition: dense_cache.hpp:48
gko::size_type
std::size_t size_type
Integral type used for allocation quantities.
Definition: types.hpp:120
gko::clone
detail::cloned_type< Pointer > clone(const Pointer &p)
Creates a unique clone of the object pointed to by p.
Definition: utils_helper.hpp:203
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
The Ginkgo namespace.
Definition: abstract_factory.hpp:48
gko::version_info::get
static const version_info & get()
Returns an instance of version_info.
Definition: version.hpp:168
gko::array< gko::stopping_status >
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::matrix::Dense::compute_norm2
void compute_norm2(ptr_param< LinOp > result) const
Computes the column-wise Euclidean (L^2) norm of this matrix.
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::write
void write(StreamType &&os, MatrixPtrType &&matrix, layout_type layout=detail::mtx_io_traits< std::remove_cv_t< detail::pointee< MatrixPtrType >>>::default_layout)
Writes a matrix into an output stream in matrix market format.
Definition: mtx_io.hpp:324
gko::log::Logger
Definition: logger.hpp:104
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::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::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::DpcppExecutor::create
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.
gko::real
constexpr auto real(const T &x)
Returns the real part of the object.
Definition: math.hpp:1041
gko::one
constexpr T one()
Returns the multiplicative identity for T.
Definition: math.hpp:803