Ginkgo Generated from branch based on main. Ginkgo version 1.11.0
A numerical linear algebra library targeting many-core architectures
Loading...
Searching...
No Matches
gko::factorization::ParIlu< ValueType, IndexType > Class Template Reference

ParILU is an incomplete LU factorization which is computed in parallel. More...

#include <ginkgo/core/factorization/par_ilu.hpp>

Inheritance diagram for gko::factorization::ParIlu< ValueType, IndexType >:
[legend]
Collaboration diagram for gko::factorization::ParIlu< ValueType, IndexType >:
[legend]

Classes

struct  parameters_type
class  Factory

Public Types

using value_type = ValueType
using index_type = IndexType
using matrix_type = matrix::Csr<ValueType, IndexType>
using l_matrix_type = matrix_type
using u_matrix_type = matrix_type
Public Types inherited from gko::Composition< default_precision >
using value_type
using transposed_type
Public Types inherited from gko::EnablePolymorphicAssignment< Composition< default_precision > >
using result_type

Public Member Functions

std::shared_ptr< const matrix_type > get_l_factor () const
std::shared_ptr< const matrix_type > get_u_factor () const
const parameters_typeget_parameters () const
Public Member Functions inherited from gko::Composition< default_precision >
const std::vector< std::shared_ptr< const LinOp > > & get_operators () const noexcept
 Returns a list of operators of the composition.
std::unique_ptr< LinOptranspose () const override
 Returns a LinOp representing the transpose of the Transposable object.
std::unique_ptr< LinOpconj_transpose () const override
 Returns a LinOp representing the conjugate transpose of the Transposable object.
Compositionoperator= (const Composition &)
 Copy-assigns a Composition.
 Composition (const Composition &)
 Copy-constructs a Composition.
Public Member Functions inherited from gko::EnableLinOp< Composition< default_precision > >
const Composition< default_precision > * apply (ptr_param< const LinOp > b, ptr_param< LinOp > x) const
Public Member Functions inherited from gko::EnablePolymorphicAssignment< Composition< default_precision > >
void convert_to (result_type *result) const override
void move_to (result_type *result) override

Static Public Member Functions

template<typename... Args>
static std::unique_ptr< Composition< ValueType > > create (Args &&... args)=delete
static auto build () -> decltype(Factory::create())
static parameters_type parse (const config::pnode &config, const config::registry &context, const config::type_descriptor &td_for_child=config::make_type_descriptor< ValueType, IndexType >())
 Create the parameters from the property_tree.
Static Public Member Functions inherited from gko::EnableCreateMethod< Composition< default_precision > >
static std::unique_ptr< Composition< default_precision > > create (Args &&... args)

Detailed Description

template<typename ValueType = default_precision, typename IndexType = int32>
class gko::factorization::ParIlu< ValueType, IndexType >

ParILU is an incomplete LU factorization which is computed in parallel.

$L$ is a lower unitriangular, while $U$ is an upper triangular matrix, which approximate a given matrix $A$ with $A \approx LU$. Here, $L$ and $U$ have the same sparsity pattern as $A$, which is also called ILU(0).

The ParILU algorithm generates the incomplete factors iteratively, using a fixed-point iteration of the form

$F(L, U) =
\begin{cases}
    \frac{1}{u_{jj}}
        \left(a_{ij}-\sum_{k=1}^{j-1}l_{ik}u_{kj}\right), \quad & i>j \\
    a_{ij}-\sum_{k=1}^{i-1}l_{ik}u_{kj}, \quad & i\leq j
\end{cases}
$

In general, the entries of $L$ and $U$ can be iterated in parallel and in asynchronous fashion, the algorithm asymptotically converges to the incomplete factors $L$ and $U$ fulfilling $\left(R = A - L \cdot
U\right)\vert_\mathcal{S} = 0\vert_\mathcal{S}$ where $\mathcal{S}$ is the pre-defined sparsity pattern (in case of ILU(0) the sparsity pattern of the system matrix $A$). The number of ParILU sweeps needed for convergence depends on the parallelism level: For sequential execution, a single sweep is sufficient, for fine-grained parallelism, the number of sweeps necessary to get a good approximation of the incomplete factors depends heavily on the problem. On the OpenMP executor, 3 sweeps usually give a decent approximation in our experiments, while GPU executors can take 10 or more iterations.

The ParILU algorithm in Ginkgo follows the design of E. Chow and A. Patel, Fine-grained Parallel Incomplete LU Factorization, SIAM Journal on Scientific Computing, 37, C169-C193 (2015).

Template Parameters
ValueTypeType of the values of all matrices used in this class
IndexTypeType of the indices of all matrices used in this class

Member Function Documentation

◆ parse()

template<typename ValueType = default_precision, typename IndexType = int32>
parameters_type gko::factorization::ParIlu< ValueType, IndexType >::parse ( const config::pnode & config,
const config::registry & context,
const config::type_descriptor & td_for_child = config::make_type_descriptor< ValueType, IndexType >() )
static

Create the parameters from the property_tree.

Because this is directly tied to the specific type, the value/index type settings within config are ignored and type_descriptor is only used for children configs.

Parameters
configthe property tree for setting
contextthe registry
td_for_childthe type descriptor for children configs. The default uses the value/index type of this class.
Returns
parameters

References gko::Composition< default_precision >::Composition().


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