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::ParIlut< ValueType, IndexType > Class Template Reference

ParILUT is an incomplete threshold-based LU factorization which is computed in parallel. More...

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

Inheritance diagram for gko::factorization::ParIlut< ValueType, IndexType >:
[legend]
Collaboration diagram for gko::factorization::ParIlut< 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::ParIlut< ValueType, IndexType >

ParILUT is an incomplete threshold-based 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 a sparsity pattern that is improved iteratively based on their element-wise magnitude. The initial sparsity pattern is chosen based on the $ILU(0)$ factorization of $A$.

One iteration of the ParILUT algorithm consists of the following steps:

  1. Calculating the residual $R = A - LU$
  2. Adding new non-zero locations from $R$ to $L$ and $U$. The new non-zero locations are initialized based on the corresponding residual value.
  3. Executing a fixed-point iteration on $L$ and $U$ according to $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}
$ For a more detailed description of the fixed-point iteration, see ParIlu.
  4. Removing the smallest entries (by magnitude) from $L$ and $U$
  5. Executing a fixed-point iteration on the (now sparser) $L$ and $U$

This ParILUT algorithm thus improves the sparsity pattern and the approximation of $L$ and $U$ simultaneously.

The implementation follows the design of H. Anzt et al., ParILUT - A Parallel Threshold ILU for GPUs, 2019 IEEE International Parallel and Distributed Processing Symposium (IPDPS), pp. 231–241.

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::ParIlut< 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: