Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
GinkgoWrappers::SolverBase< ValueType, IndexType > Class Template Reference

#include <deal.II/lac/ginkgo_solver.h>

Inheritance diagram for GinkgoWrappers::SolverBase< ValueType, IndexType >:
[legend]

Public Member Functions

 SolverBase (SolverControl &solver_control, const std::string &exec_type)
 
virtual ~SolverBase ()=default
 
void initialize (const SparseMatrix< ValueType > &matrix)
 
void apply (Vector< ValueType > &solution, const Vector< ValueType > &rhs)
 
void solve (const SparseMatrix< ValueType > &matrix, Vector< ValueType > &solution, const Vector< ValueType > &rhs)
 
SolverControlcontrol () const
 

Protected Attributes

SolverControlsolver_control
 
std::shared_ptr< gko::LinOpFactory > solver_gen
 
std::shared_ptr< gko::stop::ResidualNormReduction<>::Factory > residual_criterion
 
std::shared_ptr< gko::log::Convergence<> > convergence_logger
 
std::shared_ptr< gko::stop::Combined::Factory > combined_factory
 
std::shared_ptr< gko::Executor > executor
 

Private Member Functions

void initialize_ginkgo_log ()
 

Private Attributes

std::shared_ptr< gko::matrix::Csr< ValueType, IndexType > > system_matrix
 
const std::string exec_type
 

Detailed Description

template<typename ValueType, typename IndexType>
class GinkgoWrappers::SolverBase< ValueType, IndexType >

This class forms the base class for all of Ginkgo's iterative solvers. The various derived classes only take the additional data that is specific to them and solve the given linear system. The entire collection of solvers that Ginkgo implements is available at documentation and manual pages.

Definition at line 50 of file ginkgo_solver.h.

Constructor & Destructor Documentation

◆ SolverBase()

template<typename ValueType , typename IndexType >
SolverBase< ValueType, IndexType >::SolverBase ( SolverControl solver_control,
const std::string &  exec_type 
)

Constructor.

The exec_type defines the paradigm where the solution is computed. It is a string and the choices are "omp" , "reference" or "cuda". The respective strings create the respective executors as given below.

Ginkgo currently supports three different executor types:

  • OmpExecutor specifies that the data should be stored and the associated operations executed on an OpenMP-supporting device (e.g. host CPU);
    auto omp = gko::create<gko::OmpExecutor>();
  • CudaExecutor specifies that the data should be stored and the operations executed on the NVIDIA GPU accelerator;
    if(gko::CudaExecutor::get_num_devices() > 0 ) {
    auto cuda = gko::create<gko::CudaExecutor>();
    }
  • ReferenceExecutor executes a non-optimized reference implementation, which can be used to debug the library.
    auto ref = gko::create<gko::ReferenceExecutor>();

The following code snippet demonstrates the using of the OpenMP executor to create a solver which would use the OpenMP paradigm to the solve the system on the CPU.

auto omp = gko::create<gko::OmpExecutor>();
using cg = gko::solver::Cg<>;
auto solver_gen =
cg::build()
.with_criteria(
gko::stop::Iteration::build().with_max_iters(20u).on(omp),
gko::stop::ResidualNormReduction<>::build()
.with_reduction_factor(1e-6)
.on(omp))
.on(omp);
auto solver = solver_gen->generate(system_matrix);
solver->apply(lend(rhs), lend(solution));
std::shared_ptr< gko::matrix::Csr< ValueType, IndexType > > system_matrix
std::shared_ptr< gko::LinOpFactory > solver_gen

The solver_control object is the same as for other deal.II iterative solvers.

Definition at line 34 of file ginkgo_solver.cc.

◆ ~SolverBase()

template<typename ValueType , typename IndexType >
virtual GinkgoWrappers::SolverBase< ValueType, IndexType >::~SolverBase ( )
virtualdefault

Destructor.

Member Function Documentation

◆ initialize()

template<typename ValueType , typename IndexType >
void SolverBase< ValueType, IndexType >::initialize ( const SparseMatrix< ValueType > &  matrix)

Initialize the matrix and copy over its data to Ginkgo's data structures.

Definition at line 218 of file ginkgo_solver.cc.

◆ apply()

template<typename ValueType , typename IndexType >
void SolverBase< ValueType, IndexType >::apply ( Vector< ValueType > &  solution,
const Vector< ValueType > &  rhs 
)

Solve the linear system Ax=b. Dependent on the information provided by derived classes one of Ginkgo's linear solvers is chosen.

Definition at line 88 of file ginkgo_solver.cc.

◆ solve()

template<typename ValueType , typename IndexType >
void SolverBase< ValueType, IndexType >::solve ( const SparseMatrix< ValueType > &  matrix,
Vector< ValueType > &  solution,
const Vector< ValueType > &  rhs 
)

Solve the linear system Ax=b. Dependent on the information provided by derived classes one of Ginkgo's linear solvers is chosen.

Definition at line 292 of file ginkgo_solver.cc.

◆ control()

template<typename ValueType , typename IndexType >
SolverControl & SolverBase< ValueType, IndexType >::control

Access to the object that controls convergence.

Definition at line 209 of file ginkgo_solver.cc.

◆ initialize_ginkgo_log()

template<typename ValueType , typename IndexType >
void SolverBase< ValueType, IndexType >::initialize_ginkgo_log
private

Initialize the Ginkgo logger object with event masks. Refer to Ginkgo's logging event masks.

Definition at line 76 of file ginkgo_solver.cc.

Member Data Documentation

◆ solver_control

template<typename ValueType , typename IndexType >
SolverControl& GinkgoWrappers::SolverBase< ValueType, IndexType >::solver_control
protected

Reference to the object that controls convergence of the iterative solvers.

Definition at line 148 of file ginkgo_solver.h.

◆ solver_gen

template<typename ValueType , typename IndexType >
std::shared_ptr<gko::LinOpFactory> GinkgoWrappers::SolverBase< ValueType, IndexType >::solver_gen
protected

The Ginkgo generated solver factory object.

Definition at line 153 of file ginkgo_solver.h.

◆ residual_criterion

template<typename ValueType , typename IndexType >
std::shared_ptr<gko::stop::ResidualNormReduction<>::Factory> GinkgoWrappers::SolverBase< ValueType, IndexType >::residual_criterion
protected

The residual criterion object that controls the reduction of the residual based on the tolerance set in the solver_control member.

Definition at line 160 of file ginkgo_solver.h.

◆ convergence_logger

template<typename ValueType , typename IndexType >
std::shared_ptr<gko::log::Convergence<> > GinkgoWrappers::SolverBase< ValueType, IndexType >::convergence_logger
protected

The Ginkgo convergence logger used to check for convergence and other solver data if needed.

Definition at line 166 of file ginkgo_solver.h.

◆ combined_factory

template<typename ValueType , typename IndexType >
std::shared_ptr<gko::stop::Combined::Factory> GinkgoWrappers::SolverBase< ValueType, IndexType >::combined_factory
protected

The Ginkgo combined factory object is used to create a combined stopping criterion to be passed to the solver.

Definition at line 172 of file ginkgo_solver.h.

◆ executor

template<typename ValueType , typename IndexType >
std::shared_ptr<gko::Executor> GinkgoWrappers::SolverBase< ValueType, IndexType >::executor
protected

The execution paradigm in Ginkgo. The choices are between gko::OmpExecutor, gko::CudaExecutor and gko::ReferenceExecutor and more details can be found in Ginkgo's documentation.

Definition at line 179 of file ginkgo_solver.h.

◆ system_matrix

template<typename ValueType , typename IndexType >
std::shared_ptr<gko::matrix::Csr<ValueType, IndexType> > GinkgoWrappers::SolverBase< ValueType, IndexType >::system_matrix
private

Ginkgo matrix data structure. First template parameter is for storing the array of the non-zeros of the matrix. The second is for the row pointers and the column indices.

Todo:
Templatize based on Matrix type.

Definition at line 198 of file ginkgo_solver.h.

◆ exec_type

template<typename ValueType , typename IndexType >
const std::string GinkgoWrappers::SolverBase< ValueType, IndexType >::exec_type
private

The execution paradigm as a string to be set by the user. The choices are between omp, cuda and reference and more details can be found in Ginkgo's documentation.

Definition at line 205 of file ginkgo_solver.h.


The documentation for this class was generated from the following files: