![]() |
Ginkgo Generated from develop branch based on develop. Ginkgo version 1.8.0
A numerical linear algebra library targeting many-core architectures
|
ELL is a matrix format where stride with explicit zeros is used such that all rows have the same number of stored elements. More...
#include <ginkgo/core/matrix/ell.hpp>
Public Types | |
using | value_type = ValueType |
using | index_type = IndexType |
using | mat_data = matrix_data<ValueType, IndexType> |
using | device_mat_data = device_matrix_data<ValueType, IndexType> |
using | absolute_type = remove_complex<Ell> |
![]() | |
using | result_type |
![]() | |
using | result_type = ResultType |
![]() | |
using | value_type = ValueType |
![]() | |
using | value_type = ValueType |
using | index_type = IndexType |
![]() | |
using | value_type = ValueType |
using | index_type = IndexType |
![]() | |
using | absolute_type = AbsoluteLinOp |
Public Member Functions | |
void | convert_to (Ell< next_precision< ValueType >, IndexType > *result) const override |
void | move_to (Ell< next_precision< ValueType >, IndexType > *result) override |
void | convert_to (Dense< ValueType > *other) const override |
void | move_to (Dense< ValueType > *other) override |
void | convert_to (Csr< ValueType, IndexType > *other) const override |
void | move_to (Csr< ValueType, IndexType > *other) override |
void | read (const mat_data &data) override |
Reads a matrix from a matrix_data structure. | |
void | read (const device_mat_data &data) override |
Reads a matrix from a device_matrix_data structure. | |
void | read (device_mat_data &&data) override |
Reads a matrix from a device_matrix_data structure. | |
void | write (mat_data &data) const override |
Writes a matrix to a matrix_data structure. | |
std::unique_ptr< Diagonal< ValueType > > | extract_diagonal () const override |
Extracts the diagonal entries of the matrix into a vector. | |
std::unique_ptr< absolute_type > | compute_absolute () const override |
Gets the AbsoluteLinOp. | |
void | compute_absolute_inplace () override |
Compute absolute inplace on each element. | |
value_type * | get_values () noexcept |
Returns the values of the matrix. | |
const value_type * | get_const_values () const noexcept |
Returns the values of the matrix. | |
index_type * | get_col_idxs () noexcept |
Returns the column indexes of the matrix. | |
const index_type * | get_const_col_idxs () const noexcept |
Returns the column indexes of the matrix. | |
size_type | get_num_stored_elements_per_row () const noexcept |
Returns the number of stored elements per row. | |
size_type | get_stride () const noexcept |
Returns the stride of the matrix. | |
size_type | get_num_stored_elements () const noexcept |
Returns the number of elements explicitly stored in the matrix. | |
value_type & | val_at (size_type row, size_type idx) noexcept |
Returns the idx -th non-zero element of the row -th row . | |
value_type | val_at (size_type row, size_type idx) const noexcept |
Returns the idx -th non-zero element of the row -th row . | |
index_type & | col_at (size_type row, size_type idx) noexcept |
Returns the idx -th column index of the row -th row . | |
index_type | col_at (size_type row, size_type idx) const noexcept |
Returns the idx -th column index of the row -th row . | |
Ell & | operator= (const Ell &) |
Copy-assigns an Ell matrix. | |
Ell & | operator= (Ell &&) |
Move-assigns an Ell matrix. | |
Ell (const Ell &) | |
Copy-constructs an Ell matrix. | |
Ell (Ell &&) | |
Move-constructs an Ell matrix. | |
![]() | |
const ConcreteLinOp * | apply (ptr_param< const LinOp > b, ptr_param< LinOp > x) const |
ConcreteLinOp * | apply (ptr_param< const LinOp > b, ptr_param< LinOp > x) |
const ConcreteLinOp * | apply (ptr_param< const LinOp > alpha, ptr_param< const LinOp > b, ptr_param< const LinOp > beta, ptr_param< LinOp > x) const |
ConcreteLinOp * | apply (ptr_param< const LinOp > alpha, ptr_param< const LinOp > b, ptr_param< const LinOp > beta, ptr_param< LinOp > x) |
![]() | |
std::unique_ptr< AbstractObject > | create_default (std::shared_ptr< const Executor > exec) const |
std::unique_ptr< AbstractObject > | create_default () const |
std::unique_ptr< AbstractObject > | clone (std::shared_ptr< const Executor > exec) const |
std::unique_ptr< AbstractObject > | clone () const |
AbstractObject * | copy_from (const PolymorphicObject *other) |
template<typename Derived > | |
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, AbstractObject > * | copy_from (std::unique_ptr< Derived > &&other) |
template<typename Derived > | |
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, AbstractObject > * | copy_from (const std::unique_ptr< Derived > &other) |
AbstractObject * | copy_from (const std::shared_ptr< const PolymorphicObject > &other) |
AbstractObject * | move_from (ptr_param< PolymorphicObject > other) |
AbstractObject * | clear () |
![]() | |
PolymorphicObject & | operator= (const PolymorphicObject &) |
std::unique_ptr< PolymorphicObject > | create_default (std::shared_ptr< const Executor > exec) const |
Creates a new "default" object of the same dynamic type as this object. | |
std::unique_ptr< PolymorphicObject > | create_default () const |
Creates a new "default" object of the same dynamic type as this object. | |
std::unique_ptr< PolymorphicObject > | clone (std::shared_ptr< const Executor > exec) const |
Creates a clone of the object. | |
std::unique_ptr< PolymorphicObject > | clone () const |
Creates a clone of the object. | |
PolymorphicObject * | copy_from (const PolymorphicObject *other) |
Copies another object into this object. | |
template<typename Derived , typename Deleter > | |
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, PolymorphicObject > * | copy_from (std::unique_ptr< Derived, Deleter > &&other) |
Moves another object into this object. | |
template<typename Derived , typename Deleter > | |
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, PolymorphicObject > * | copy_from (const std::unique_ptr< Derived, Deleter > &other) |
Copies another object into this object. | |
PolymorphicObject * | copy_from (const std::shared_ptr< const PolymorphicObject > &other) |
Copies another object into this object. | |
PolymorphicObject * | move_from (ptr_param< PolymorphicObject > other) |
Moves another object into this object. | |
PolymorphicObject * | clear () |
Transforms the object into its default state. | |
std::shared_ptr< const Executor > | get_executor () const noexcept |
Returns the Executor of the object. | |
![]() | |
void | add_logger (std::shared_ptr< const Logger > logger) override |
Adds a new logger to the list of subscribed loggers. | |
void | remove_logger (const Logger *logger) override |
Removes a logger from the list of subscribed loggers. | |
void | remove_logger (ptr_param< const Logger > logger) |
const std::vector< std::shared_ptr< const Logger > > & | get_loggers () const override |
Returns the vector containing all loggers registered at this object. | |
void | clear_loggers () override |
Remove all loggers registered at this object. | |
![]() | |
void | remove_logger (ptr_param< const Logger > logger) |
![]() | |
void | convert_to (result_type *result) const override |
Converts the implementer to an object of type result_type. | |
void | move_to (result_type *result) override |
Converts the implementer to an object of type result_type by moving data from this object. | |
![]() | |
virtual void | convert_to (result_type *result) const =0 |
Converts the implementer to an object of type result_type. | |
void | convert_to (ptr_param< result_type > result) const |
virtual void | move_to (result_type *result)=0 |
Converts the implementer to an object of type result_type by moving data from this object. | |
void | move_to (ptr_param< result_type > result) |
![]() | |
std::unique_ptr< LinOp > | extract_diagonal_linop () const override |
Extracts the diagonal entries of the matrix into a vector. | |
![]() | |
void | read (const matrix_assembly_data< ValueType, IndexType > &data) |
Reads a matrix from a matrix_assembly_data structure. | |
![]() | |
std::unique_ptr< LinOp > | compute_absolute_linop () const override |
Gets the absolute LinOp. | |
Static Public Member Functions | |
static std::unique_ptr< Ell > | create (std::shared_ptr< const Executor > exec, const dim< 2 > &size={}, size_type num_stored_elements_per_row=0, size_type stride=0) |
Creates an uninitialized Ell matrix of the specified size. | |
static std::unique_ptr< Ell > | create (std::shared_ptr< const Executor > exec, const dim< 2 > &size, array< value_type > values, array< index_type > col_idxs, size_type num_stored_elements_per_row, size_type stride) |
Creates an ELL matrix from already allocated (and initialized) column index and value arrays. | |
template<typename InputValueType , typename InputColumnIndexType > | |
static std::unique_ptr< Ell > | create (std::shared_ptr< const Executor > exec, const dim< 2 > &size, std::initializer_list< InputValueType > values, std::initializer_list< InputColumnIndexType > col_idxs, size_type num_stored_elements_per_row, size_type stride) |
create(std::shared_ptr<const Executor>,const dim<2>&, array<value_type>, array<index_type>, size_type,size_type) | |
static std::unique_ptr< const Ell > | create_const (std::shared_ptr< const Executor > exec, const dim< 2 > &size, gko::detail::const_array_view< ValueType > &&values, gko::detail::const_array_view< IndexType > &&col_idxs, size_type num_stored_elements_per_row, size_type stride) |
Creates a constant (immutable) Ell matrix from a set of constant arrays. | |
Friends | |
class | EnablePolymorphicObject< Ell, LinOp > |
class | Dense< ValueType > |
class | Coo< ValueType, IndexType > |
class | Csr< ValueType, IndexType > |
class | Ell< to_complex< ValueType >, IndexType > |
class | Ell< next_precision< ValueType >, IndexType > |
class | Hybrid< ValueType, IndexType > |
ELL is a matrix format where stride with explicit zeros is used such that all rows have the same number of stored elements.
The number of elements stored in each row is the largest number of nonzero elements in any of the rows (obtainable through get_num_stored_elements_per_row() method). This removes the need of a row pointer like in the CSR format, and allows for SIMD processing of the distinct rows. For efficient processing, the nonzero elements and the corresponding column indices are stored in column-major fashion. The columns are padded to the length by user-defined stride parameter whose default value is the number of rows of the matrix.
This implementation uses the column index value invalid_index<IndexType>() to mark padding entries that are not part of the sparsity pattern.
ValueType | precision of matrix elements |
IndexType | precision of matrix indexes |
gko::matrix::Ell< ValueType, IndexType >::Ell | ( | const Ell< ValueType, IndexType > & | ) |
Copy-constructs an Ell matrix.
Inherits executor and dimensions, but copies data without padding.
gko::matrix::Ell< ValueType, IndexType >::Ell | ( | Ell< ValueType, IndexType > && | ) |
Move-constructs an Ell matrix.
Inherits executor, dimensions and data with padding. The moved-from object is empty (0x0 with empty Array).
|
inlinenoexcept |
Returns the idx
-th column index of the row
-th row .
row | the row of the requested element |
idx | the idx-th stored element of the row |
References gko::matrix::Ell< ValueType, IndexType >::get_const_col_idxs(), and gko::one().
|
inlinenoexcept |
Returns the idx
-th column index of the row
-th row .
row | the row of the requested element |
idx | the idx-th stored element of the row |
References gko::matrix::Ell< ValueType, IndexType >::get_col_idxs(), and gko::one().
|
overridevirtual |
Gets the AbsoluteLinOp.
Implements gko::EnableAbsoluteComputation< AbsoluteLinOp >.
|
overridevirtual |
Compute absolute inplace on each element.
Implements gko::AbsoluteComputable.
|
static |
Creates an ELL matrix from already allocated (and initialized) column index and value arrays.
exec | Executor associated to the matrix |
size | size of the matrix |
values | array of matrix values |
col_idxs | array of column indexes |
num_stored_elements_per_row | the number of stored elements per row |
stride | stride of the rows |
col_idxs
or values
is not an rvalue, not an array of IndexType and ValueType, respectively, or is on the wrong executor, an internal copy of that array will be created, and the original array data will not be used in the matrix.
|
inlinestatic |
References gko::matrix::Ell< ValueType, IndexType >::create(), and gko::one().
|
static |
Creates an uninitialized Ell matrix of the specified size.
exec | Executor associated to the matrix |
size | size of the matrix |
num_stored_elements_per_row | the number of stored elements per row |
stride | stride of the columns. If it is set to 0, size[0] will be used instead. |
Referenced by gko::matrix::Ell< ValueType, IndexType >::create().
|
static |
Creates a constant (immutable) Ell matrix from a set of constant arrays.
exec | the executor to create the matrix on |
size | the dimensions of the matrix |
values | the value array of the matrix |
col_idxs | the column index array of the matrix |
num_stored_elements_per_row | the number of stored nonzeros per row |
stride | the column-stride of the value and column index array |
|
overridevirtual |
Extracts the diagonal entries of the matrix into a vector.
diag | the vector into which the diagonal will be written |
Implements gko::DiagonalExtractable< ValueType >.
|
inlinenoexcept |
Returns the column indexes of the matrix.
References gko::array< ValueType >::get_data().
Referenced by gko::matrix::Ell< ValueType, IndexType >::col_at().
|
inlinenoexcept |
Returns the column indexes of the matrix.
References gko::array< ValueType >::get_const_data().
Referenced by gko::matrix::Ell< ValueType, IndexType >::col_at().
|
inlinenoexcept |
Returns the values of the matrix.
References gko::array< ValueType >::get_const_data().
|
inlinenoexcept |
Returns the number of elements explicitly stored in the matrix.
References gko::array< ValueType >::get_size().
|
inlinenoexcept |
Returns the number of stored elements per row.
|
inlinenoexcept |
Returns the stride of the matrix.
|
inlinenoexcept |
Returns the values of the matrix.
References gko::array< ValueType >::get_data().
Ell & gko::matrix::Ell< ValueType, IndexType >::operator= | ( | const Ell< ValueType, IndexType > & | ) |
Copy-assigns an Ell matrix.
Preserves the executor, reallocates the matrix with minimal stride if the dimensions don't match, then copies the data over, ignoring padding.
Ell & gko::matrix::Ell< ValueType, IndexType >::operator= | ( | Ell< ValueType, IndexType > && | ) |
Move-assigns an Ell matrix.
Preserves the executor, moves the data over preserving size and stride. Leaves the moved-from object in an empty state (0x0 with empty Array).
|
overridevirtual |
Reads a matrix from a device_matrix_data structure.
data | the device_matrix_data structure. |
Reimplemented from gko::ReadableFromMatrixData< ValueType, IndexType >.
|
overridevirtual |
Reads a matrix from a matrix_data structure.
data | the matrix_data structure |
Implements gko::ReadableFromMatrixData< ValueType, IndexType >.
|
overridevirtual |
Reads a matrix from a device_matrix_data structure.
The structure may be emptied by this function.
data | the device_matrix_data structure. |
Reimplemented from gko::ReadableFromMatrixData< ValueType, IndexType >.
|
inlinenoexcept |
Returns the idx
-th non-zero element of the row
-th row .
row | the row of the requested element |
idx | the idx-th stored element of the row |
References gko::array< ValueType >::get_const_data(), and gko::one().
|
inlinenoexcept |
Returns the idx
-th non-zero element of the row
-th row .
row | the row of the requested element |
idx | the idx-th stored element of the row |
References gko::array< ValueType >::get_data(), and gko::one().
|
overridevirtual |
Writes a matrix to a matrix_data structure.
data | the matrix_data structure |
Implements gko::WritableToMatrixData< ValueType, IndexType >.