![]() |
Ginkgo
Generated from pipelines/1744748943 branch based on develop. Ginkgo version 1.10.0
A numerical linear algebra library targeting many-core architectures
|
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>
Classes | |
class | Factory |
struct | parameters_type |
Public Types | |
using | value_type = ValueType |
using | index_type = IndexType |
using | matrix_type = batch::matrix::Csr< ValueType, IndexType > |
![]() | |
using | result_type = 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_type & | get_parameters () const |
![]() | |
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 () |
![]() | |
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... | |
![]() | |
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 > |
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.
ValueType | value precision of matrix elements |
IndexType | index precision of matrix elements |
|
inlinenoexcept |
Returns the block pointers.
References gko::array< ValueType >::get_const_data().
|
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
References gko::array< ValueType >::get_const_data().
|
inlinenoexcept |
Returns the cumulative blocks storage array.
References gko::array< ValueType >::get_const_data().
|
inlinenoexcept |
Returns the mapping between the blocks and the row id.
References gko::array< ValueType >::get_const_data().
|
inlinenoexcept |
Returns the max block size.
References gko::batch::preconditioner::Jacobi< ValueType, IndexType >::parameters_type::max_block_size.
|
inlinenoexcept |
Returns the number of blocks in an individual batch entry.
|
inlinenoexcept |
Returns the number of elements explicitly stored in the dense blocks.
References gko::array< ValueType >::get_size().