Ginkgo  Generated from pipelines/1330831941 branch based on master. Ginkgo version 1.8.0
A numerical linear algebra library targeting many-core architectures
Classes | Public Types | Public Member Functions | Static Public Member Functions | Friends | List of all members
gko::batch::preconditioner::Jacobi< ValueType, IndexType > Class Template Referencefinal

A block-Jacobi preconditioner is a block-diagonal linear operator, obtained by inverting the diagonal blocks (stored in a dense row major fashion) of the source operator. More...

#include <ginkgo/core/preconditioner/batch_jacobi.hpp>

Inheritance diagram for gko::batch::preconditioner::Jacobi< ValueType, IndexType >:
[legend]
Collaboration diagram for gko::batch::preconditioner::Jacobi< ValueType, IndexType >:
[legend]

Classes

class  Factory
 
struct  parameters_type
 

Public Types

using value_type = ValueType
 
using index_type = IndexType
 
using matrix_type = batch::matrix::Csr< ValueType, IndexType >
 
- Public Types inherited from gko::EnablePolymorphicAssignment< Jacobi< ValueType, IndexType > >
using result_type = Jacobi< ValueType, IndexType >
 
- Public Types inherited from gko::ConvertibleTo< Jacobi< ValueType, IndexType > >
using result_type = Jacobi< ValueType, IndexType >
 

Public Member Functions

const index_type * get_const_block_pointers () const noexcept
 Returns the block pointers. More...
 
const index_type * get_const_map_block_to_row () const noexcept
 Returns the mapping between the blocks and the row id. More...
 
const index_type * get_const_blocks_cumulative_offsets () const noexcept
 Returns the cumulative blocks storage array. More...
 
uint32 get_max_block_size () const noexcept
 Returns the max block size. More...
 
size_type get_num_blocks () const noexcept
 Returns the number of blocks in an individual batch entry. More...
 
const value_type * get_const_blocks () const noexcept
 Returns the pointer to the memory used for storing the block data. More...
 
size_type get_num_stored_elements () const noexcept
 Returns the number of elements explicitly stored in the dense blocks. More...
 
const parameters_typeget_parameters () const
 
- Public Member Functions inherited from gko::EnableAbstractPolymorphicObject< Jacobi< ValueType, IndexType >, BatchLinOp >
std::unique_ptr< Jacobi< ValueType, IndexType > > create_default (std::shared_ptr< const Executor > exec) const
 
std::unique_ptr< Jacobi< ValueType, IndexType > > create_default () const
 
std::unique_ptr< Jacobi< ValueType, IndexType > > clone (std::shared_ptr< const Executor > exec) const
 
std::unique_ptr< Jacobi< ValueType, IndexType > > clone () const
 
Jacobi< ValueType, IndexType > * copy_from (const PolymorphicObject *other)
 
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, Jacobi< ValueType, IndexType > > * copy_from (std::unique_ptr< Derived > &&other)
 
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, Jacobi< ValueType, IndexType > > * copy_from (const std::unique_ptr< Derived > &other)
 
Jacobi< ValueType, IndexType > * copy_from (const std::shared_ptr< const PolymorphicObject > &other)
 
Jacobi< ValueType, IndexType > * move_from (ptr_param< PolymorphicObject > other)
 
Jacobi< ValueType, IndexType > * clear ()
 
- Public Member Functions inherited from gko::EnablePolymorphicAssignment< Jacobi< ValueType, IndexType > >
void convert_to (result_type *result) const override
 Converts the implementer to an object of type result_type. More...
 
void move_to (result_type *result) override
 Converts the implementer to an object of type result_type by moving data from this object. More...
 
- Public Member Functions inherited from gko::ConvertibleTo< Jacobi< ValueType, IndexType > >
void convert_to (ptr_param< result_type > result) const
 
void move_to (ptr_param< result_type > result)
 

Static Public Member Functions

static auto build () -> decltype(Factory ::create())
 

Friends

class EnableBatchLinOp< Jacobi >
 
class EnablePolymorphicObject< Jacobi, BatchLinOp >
 

Detailed Description

template<typename ValueType = default_precision, typename IndexType = int32>
class gko::batch::preconditioner::Jacobi< ValueType, IndexType >

A block-Jacobi preconditioner is a block-diagonal linear operator, obtained by inverting the diagonal blocks (stored in a dense row major fashion) of the source operator.

With the batched preconditioners, it is required that all items in the batch have the same sparsity pattern. The detection of the blocks and the block pointers require that the sparsity pattern of all the items be the same. Other cases is undefined behaviour. The input batch matrix must be in batch::Csr matrix format or must be convertible to batch::Csr matrix format. The block detection algorithm and the conversion to dense blocks kernels require this assumption.

Note
In a fashion similar to the non-batched Jacobi preconditioner, the maximum possible size of the diagonal blocks is equal to the maximum warp size on the device (32 for NVIDIA GPUs, 64 for AMD GPUs).
Template Parameters
ValueTypevalue precision of matrix elements
IndexTypeindex precision of matrix elements

Member Function Documentation

◆ get_const_block_pointers()

template<typename ValueType = default_precision, typename IndexType = int32>
const index_type* gko::batch::preconditioner::Jacobi< ValueType, IndexType >::get_const_block_pointers ( ) const
inlinenoexcept

Returns the block pointers.

Note
Returns nullptr in case of a scalar jacobi preconditioner (max_block_size = 1).
Returns
the block pointers

References gko::array< ValueType >::get_const_data().

◆ get_const_blocks()

template<typename ValueType = default_precision, typename IndexType = int32>
const value_type* gko::batch::preconditioner::Jacobi< ValueType, IndexType >::get_const_blocks ( ) const
inlinenoexcept

Returns the pointer to the memory used for storing the block data.

Element (i, j) of the block, which belongs to the batch entry with index = "batch_id" and has local id = "block_id" within its batch entry is stored at the address = get_const_blocks() + detail::get_global_block_offset(batch_id, num_blocks, block_id, cumulative_blocks_storage) + i * detail::get_stride(block_id, block_pointers) + j

Note
Returns nullptr in case of a scalar jacobi preconditioner (max_block_size = 1). The blocks array is empty in case of scalar jacobi preconditioner as the preconditioner is generated inside the batched solver kernel.
Returns
the pointer to the memory used for storing the block data

References gko::array< ValueType >::get_const_data().

◆ get_const_blocks_cumulative_offsets()

template<typename ValueType = default_precision, typename IndexType = int32>
const index_type* gko::batch::preconditioner::Jacobi< ValueType, IndexType >::get_const_blocks_cumulative_offsets ( ) const
inlinenoexcept

Returns the cumulative blocks storage array.

Note
Returns nullptr in case of a scalar jacobi preconditioner (max_block_size = 1).

References gko::array< ValueType >::get_const_data().

◆ get_const_map_block_to_row()

template<typename ValueType = default_precision, typename IndexType = int32>
const index_type* gko::batch::preconditioner::Jacobi< ValueType, IndexType >::get_const_map_block_to_row ( ) const
inlinenoexcept

Returns the mapping between the blocks and the row id.

Note
Returns nullptr in case of a scalar jacobi preconditioner (max_block_size = 1).

References gko::array< ValueType >::get_const_data().

◆ get_max_block_size()

template<typename ValueType = default_precision, typename IndexType = int32>
uint32 gko::batch::preconditioner::Jacobi< ValueType, IndexType >::get_max_block_size ( ) const
inlinenoexcept

Returns the max block size.

Returns
the max block size

References gko::batch::preconditioner::Jacobi< ValueType, IndexType >::parameters_type::max_block_size.

◆ get_num_blocks()

template<typename ValueType = default_precision, typename IndexType = int32>
size_type gko::batch::preconditioner::Jacobi< ValueType, IndexType >::get_num_blocks ( ) const
inlinenoexcept

Returns the number of blocks in an individual batch entry.

Returns
the number of blocks in an individual batch entry.

◆ get_num_stored_elements()

template<typename ValueType = default_precision, typename IndexType = int32>
size_type gko::batch::preconditioner::Jacobi< ValueType, IndexType >::get_num_stored_elements ( ) const
inlinenoexcept

Returns the number of elements explicitly stored in the dense blocks.

Note
Returns 0 in case of scalar jacobi preconditioner as the preconditioner is generated inside the batched solver kernels, hence, blocks array storage is not required.
Returns
the number of elements explicitly stored in the dense blocks.

References gko::array< ValueType >::get_size().


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