Ginkgo
Generated from pipelines/1330831941 branch based on master. Ginkgo version 1.8.0
A numerical linear algebra library targeting many-core architectures
|
CSR is a matrix format which stores only the nonzero coefficients by compressing each row of the matrix (compressed sparse row format). More...
#include <ginkgo/core/matrix/csr.hpp>
Classes | |
class | automatical |
class | classical |
classical is a strategy_type which uses the same number of threads on each row. More... | |
class | cusparse |
cusparse is a strategy_type which uses the sparselib csr. More... | |
class | load_balance |
load_balance is a strategy_type which uses the load balance algorithm. More... | |
class | merge_path |
merge_path is a strategy_type which uses the merge_path algorithm. More... | |
class | sparselib |
sparselib is a strategy_type which uses the sparselib csr. More... | |
class | strategy_type |
strategy_type is to decide how to set the csr algorithm. More... | |
Public Types | |
using | value_type = ValueType |
using | index_type = IndexType |
using | transposed_type = Csr< ValueType, IndexType > |
using | mat_data = matrix_data< ValueType, IndexType > |
using | device_mat_data = device_matrix_data< ValueType, IndexType > |
using | absolute_type = remove_complex< Csr > |
Public Types inherited from gko::EnablePolymorphicAssignment< Csr< ValueType, IndexType > > | |
using | result_type = Csr< ValueType, IndexType > |
Public Types inherited from gko::ConvertibleTo< Csr< ValueType, IndexType > > | |
using | result_type = Csr< ValueType, IndexType > |
Public Types inherited from gko::ConvertibleTo< Csr< next_precision< ValueType >, IndexType > > | |
using | result_type = Csr< next_precision< ValueType >, IndexType > |
Public Types inherited from gko::ConvertibleTo< Dense< ValueType > > | |
using | result_type = Dense< ValueType > |
Public Types inherited from gko::ConvertibleTo< Coo< ValueType, IndexType > > | |
using | result_type = Coo< ValueType, IndexType > |
Public Types inherited from gko::ConvertibleTo< Ell< ValueType, IndexType > > | |
using | result_type = Ell< ValueType, IndexType > |
Public Types inherited from gko::ConvertibleTo< Fbcsr< ValueType, IndexType > > | |
using | result_type = Fbcsr< ValueType, IndexType > |
Public Types inherited from gko::ConvertibleTo< Hybrid< ValueType, IndexType > > | |
using | result_type = Hybrid< ValueType, IndexType > |
Public Types inherited from gko::ConvertibleTo< Sellp< ValueType, IndexType > > | |
using | result_type = Sellp< ValueType, IndexType > |
Public Types inherited from gko::ConvertibleTo< SparsityCsr< ValueType, IndexType > > | |
using | result_type = SparsityCsr< ValueType, IndexType > |
Public Types inherited from gko::DiagonalExtractable< ValueType > | |
using | value_type = ValueType |
Public Types inherited from gko::ReadableFromMatrixData< ValueType, IndexType > | |
using | value_type = ValueType |
using | index_type = IndexType |
Public Types inherited from gko::WritableToMatrixData< ValueType, IndexType > | |
using | value_type = ValueType |
using | index_type = IndexType |
Public Types inherited from gko::EnableAbsoluteComputation< remove_complex< Csr< ValueType, IndexType > > > | |
using | absolute_type = remove_complex< Csr< ValueType, IndexType > > |
Public Member Functions | |
void | convert_to (Csr< next_precision< ValueType >, IndexType > *result) const override |
void | move_to (Csr< next_precision< ValueType >, IndexType > *result) override |
void | convert_to (Dense< ValueType > *other) const override |
void | move_to (Dense< ValueType > *other) override |
void | convert_to (Coo< ValueType, IndexType > *result) const override |
void | move_to (Coo< ValueType, IndexType > *result) override |
void | convert_to (Ell< ValueType, IndexType > *result) const override |
void | move_to (Ell< ValueType, IndexType > *result) override |
void | convert_to (Fbcsr< ValueType, IndexType > *result) const override |
void | move_to (Fbcsr< ValueType, IndexType > *result) override |
void | convert_to (Hybrid< ValueType, IndexType > *result) const override |
void | move_to (Hybrid< ValueType, IndexType > *result) override |
void | convert_to (Sellp< ValueType, IndexType > *result) const override |
void | move_to (Sellp< ValueType, IndexType > *result) override |
void | convert_to (SparsityCsr< ValueType, IndexType > *result) const override |
void | move_to (SparsityCsr< ValueType, IndexType > *result) override |
void | read (const mat_data &data) override |
Reads a matrix from a matrix_data structure. More... | |
void | read (const device_mat_data &data) override |
Reads a matrix from a device_matrix_data structure. More... | |
void | read (device_mat_data &&data) override |
Reads a matrix from a device_matrix_data structure. More... | |
void | write (mat_data &data) const override |
Writes a matrix to a matrix_data structure. More... | |
std::unique_ptr< LinOp > | transpose () const override |
Returns a LinOp representing the transpose of the Transposable object. More... | |
std::unique_ptr< LinOp > | conj_transpose () const override |
Returns a LinOp representing the conjugate transpose of the Transposable object. More... | |
std::unique_ptr< Csr > | permute (ptr_param< const Permutation< index_type >> permutation, permute_mode mode=permute_mode::symmetric) const |
Creates a permuted copy of this matrix with the given permutation . More... | |
std::unique_ptr< Csr > | permute (ptr_param< const Permutation< index_type >> row_permutation, ptr_param< const Permutation< index_type >> column_permutation, bool invert=false) const |
Creates a non-symmetrically permuted copy of this matrix with the given row and column permutations and . More... | |
std::unique_ptr< Csr > | scale_permute (ptr_param< const ScaledPermutation< value_type, index_type >> permutation, permute_mode=permute_mode::symmetric) const |
Creates a scaled and permuted copy of this matrix. More... | |
std::unique_ptr< Csr > | scale_permute (ptr_param< const ScaledPermutation< value_type, index_type >> row_permutation, ptr_param< const ScaledPermutation< value_type, index_type >> column_permutation, bool invert=false) const |
Creates a scaled and permuted copy of this matrix. More... | |
std::unique_ptr< LinOp > | permute (const array< IndexType > *permutation_indices) const override |
Returns a LinOp representing the symmetric row and column permutation of the Permutable object. More... | |
std::unique_ptr< LinOp > | inverse_permute (const array< IndexType > *inverse_permutation_indices) const override |
Returns a LinOp representing the symmetric inverse row and column permutation of the Permutable object. More... | |
std::unique_ptr< LinOp > | row_permute (const array< IndexType > *permutation_indices) const override |
Returns a LinOp representing the row permutation of the Permutable object. More... | |
std::unique_ptr< LinOp > | column_permute (const array< IndexType > *permutation_indices) const override |
Returns a LinOp representing the column permutation of the Permutable object. More... | |
std::unique_ptr< LinOp > | inverse_row_permute (const array< IndexType > *inverse_permutation_indices) const override |
Returns a LinOp representing the row permutation of the inverse permuted object. More... | |
std::unique_ptr< LinOp > | inverse_column_permute (const array< IndexType > *inverse_permutation_indices) const override |
Returns a LinOp representing the row permutation of the inverse permuted object. More... | |
std::unique_ptr< Diagonal< ValueType > > | extract_diagonal () const override |
Extracts the diagonal entries of the matrix into a vector. More... | |
std::unique_ptr< absolute_type > | compute_absolute () const override |
Gets the AbsoluteLinOp. More... | |
void | compute_absolute_inplace () override |
Compute absolute inplace on each element. | |
void | sort_by_column_index () |
Sorts all (value, col_idx) pairs in each row by column index. | |
bool | is_sorted_by_column_index () const |
value_type * | get_values () noexcept |
Returns the values of the matrix. More... | |
const value_type * | get_const_values () const noexcept |
Returns the values of the matrix. More... | |
index_type * | get_col_idxs () noexcept |
Returns the column indexes of the matrix. More... | |
const index_type * | get_const_col_idxs () const noexcept |
Returns the column indexes of the matrix. More... | |
index_type * | get_row_ptrs () noexcept |
Returns the row pointers of the matrix. More... | |
const index_type * | get_const_row_ptrs () const noexcept |
Returns the row pointers of the matrix. More... | |
index_type * | get_srow () noexcept |
Returns the starting rows. More... | |
const index_type * | get_const_srow () const noexcept |
Returns the starting rows. More... | |
size_type | get_num_srow_elements () const noexcept |
Returns the number of the srow stored elements (involved warps) More... | |
size_type | get_num_stored_elements () const noexcept |
Returns the number of elements explicitly stored in the matrix. More... | |
std::shared_ptr< strategy_type > | get_strategy () const noexcept |
Returns the strategy. More... | |
void | set_strategy (std::shared_ptr< strategy_type > strategy) |
Set the strategy. More... | |
void | scale (ptr_param< const LinOp > alpha) |
Scales the matrix with a scalar. More... | |
void | inv_scale (ptr_param< const LinOp > alpha) |
Scales the matrix with the inverse of a scalar. More... | |
std::unique_ptr< Csr< ValueType, IndexType > > | create_submatrix (const index_set< IndexType > &row_index_set, const index_set< IndexType > &column_index_set) const |
Creates a submatrix from this Csr matrix given row and column index_set objects. More... | |
std::unique_ptr< Csr< ValueType, IndexType > > | create_submatrix (const span &row_span, const span &column_span) const |
Creates a submatrix from this Csr matrix given row and column spans. More... | |
Csr & | operator= (const Csr &) |
Copy-assigns a Csr matrix. More... | |
Csr & | operator= (Csr &&) |
Move-assigns a Csr matrix. More... | |
Csr (const Csr &) | |
Copy-constructs a Csr matrix. More... | |
Csr (Csr &&) | |
Move-constructs a Csr matrix. More... | |
Public Member Functions inherited from gko::EnableLinOp< Csr< ValueType, IndexType > > | |
const Csr< ValueType, IndexType > * | apply (ptr_param< const LinOp > b, ptr_param< LinOp > x) const |
Csr< ValueType, IndexType > * | apply (ptr_param< const LinOp > b, ptr_param< LinOp > x) |
const Csr< ValueType, IndexType > * | apply (ptr_param< const LinOp > alpha, ptr_param< const LinOp > b, ptr_param< const LinOp > beta, ptr_param< LinOp > x) const |
Csr< ValueType, IndexType > * | apply (ptr_param< const LinOp > alpha, ptr_param< const LinOp > b, ptr_param< const LinOp > beta, ptr_param< LinOp > x) |
Public Member Functions inherited from gko::EnableAbstractPolymorphicObject< Csr< ValueType, IndexType >, LinOp > | |
std::unique_ptr< Csr< ValueType, IndexType > > | create_default (std::shared_ptr< const Executor > exec) const |
std::unique_ptr< Csr< ValueType, IndexType > > | create_default () const |
std::unique_ptr< Csr< ValueType, IndexType > > | clone (std::shared_ptr< const Executor > exec) const |
std::unique_ptr< Csr< ValueType, IndexType > > | clone () const |
Csr< ValueType, IndexType > * | copy_from (const PolymorphicObject *other) |
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, Csr< ValueType, IndexType > > * | copy_from (std::unique_ptr< Derived > &&other) |
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, Csr< ValueType, IndexType > > * | copy_from (const std::unique_ptr< Derived > &other) |
Csr< ValueType, IndexType > * | copy_from (const std::shared_ptr< const PolymorphicObject > &other) |
Csr< ValueType, IndexType > * | move_from (ptr_param< PolymorphicObject > other) |
Csr< ValueType, IndexType > * | clear () |
Public Member Functions inherited from gko::LinOp | |
LinOp * | apply (ptr_param< const LinOp > b, ptr_param< LinOp > x) |
Applies a linear operator to a vector (or a sequence of vectors). More... | |
const LinOp * | apply (ptr_param< const LinOp > b, ptr_param< LinOp > x) const |
LinOp * | apply (ptr_param< const LinOp > alpha, ptr_param< const LinOp > b, ptr_param< const LinOp > beta, ptr_param< LinOp > x) |
Performs the operation x = alpha * op(b) + beta * x. More... | |
const LinOp * | apply (ptr_param< const LinOp > alpha, ptr_param< const LinOp > b, ptr_param< const LinOp > beta, ptr_param< LinOp > x) const |
const dim< 2 > & | get_size () const noexcept |
Returns the size of the operator. More... | |
virtual bool | apply_uses_initial_guess () const |
Returns true if the linear operator uses the data given in x as an initial guess. More... | |
LinOp & | operator= (const LinOp &)=default |
Copy-assigns a LinOp. More... | |
LinOp & | operator= (LinOp &&other) |
Move-assigns a LinOp. More... | |
LinOp (const LinOp &)=default | |
Copy-constructs a LinOp. More... | |
LinOp (LinOp &&other) | |
Move-constructs a LinOp. More... | |
Public Member Functions inherited from gko::EnableAbstractPolymorphicObject< LinOp > | |
std::unique_ptr< LinOp > | create_default (std::shared_ptr< const Executor > exec) const |
std::unique_ptr< LinOp > | create_default () const |
std::unique_ptr< LinOp > | clone (std::shared_ptr< const Executor > exec) const |
std::unique_ptr< LinOp > | clone () const |
LinOp * | copy_from (const PolymorphicObject *other) |
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, LinOp > * | copy_from (std::unique_ptr< Derived > &&other) |
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, LinOp > * | copy_from (const std::unique_ptr< Derived > &other) |
LinOp * | copy_from (const std::shared_ptr< const PolymorphicObject > &other) |
LinOp * | move_from (ptr_param< PolymorphicObject > other) |
LinOp * | clear () |
Public Member Functions inherited from gko::PolymorphicObject | |
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. More... | |
std::unique_ptr< PolymorphicObject > | create_default () const |
Creates a new "default" object of the same dynamic type as this object. More... | |
std::unique_ptr< PolymorphicObject > | clone (std::shared_ptr< const Executor > exec) const |
Creates a clone of the object. More... | |
std::unique_ptr< PolymorphicObject > | clone () const |
Creates a clone of the object. More... | |
PolymorphicObject * | copy_from (const PolymorphicObject *other) |
Copies another object into this object. More... | |
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. More... | |
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. More... | |
PolymorphicObject * | copy_from (const std::shared_ptr< const PolymorphicObject > &other) |
Copies another object into this object. More... | |
PolymorphicObject * | move_from (ptr_param< PolymorphicObject > other) |
Moves another object into this object. More... | |
PolymorphicObject * | clear () |
Transforms the object into its default state. More... | |
std::shared_ptr< const Executor > | get_executor () const noexcept |
Returns the Executor of the object. More... | |
Public Member Functions inherited from gko::log::EnableLogging< PolymorphicObject > | |
void | add_logger (std::shared_ptr< const Logger > logger) override |
void | remove_logger (const Logger *logger) override |
void | remove_logger (ptr_param< const Logger > logger) |
const std::vector< std::shared_ptr< const Logger > > & | get_loggers () const override |
void | clear_loggers () override |
Public Member Functions inherited from gko::log::Loggable | |
void | remove_logger (ptr_param< const Logger > logger) |
Public Member Functions inherited from gko::EnablePolymorphicAssignment< Csr< 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< Csr< ValueType, IndexType > > | |
void | convert_to (ptr_param< result_type > result) const |
void | move_to (ptr_param< result_type > result) |
Public Member Functions inherited from gko::ConvertibleTo< Csr< next_precision< ValueType >, IndexType > > | |
virtual void | convert_to (result_type *result) const=0 |
Converts the implementer to an object of type result_type. More... | |
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. More... | |
void | move_to (ptr_param< result_type > result) |
Public Member Functions inherited from gko::ConvertibleTo< Dense< ValueType > > | |
virtual void | convert_to (result_type *result) const=0 |
Converts the implementer to an object of type result_type. More... | |
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. More... | |
void | move_to (ptr_param< result_type > result) |
Public Member Functions inherited from gko::ConvertibleTo< Coo< ValueType, IndexType > > | |
virtual void | convert_to (result_type *result) const=0 |
Converts the implementer to an object of type result_type. More... | |
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. More... | |
void | move_to (ptr_param< result_type > result) |
Public Member Functions inherited from gko::ConvertibleTo< Ell< ValueType, IndexType > > | |
virtual void | convert_to (result_type *result) const=0 |
Converts the implementer to an object of type result_type. More... | |
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. More... | |
void | move_to (ptr_param< result_type > result) |
Public Member Functions inherited from gko::ConvertibleTo< Fbcsr< ValueType, IndexType > > | |
virtual void | convert_to (result_type *result) const=0 |
Converts the implementer to an object of type result_type. More... | |
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. More... | |
void | move_to (ptr_param< result_type > result) |
Public Member Functions inherited from gko::ConvertibleTo< Hybrid< ValueType, IndexType > > | |
virtual void | convert_to (result_type *result) const=0 |
Converts the implementer to an object of type result_type. More... | |
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. More... | |
void | move_to (ptr_param< result_type > result) |
Public Member Functions inherited from gko::ConvertibleTo< Sellp< ValueType, IndexType > > | |
virtual void | convert_to (result_type *result) const=0 |
Converts the implementer to an object of type result_type. More... | |
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. More... | |
void | move_to (ptr_param< result_type > result) |
Public Member Functions inherited from gko::ConvertibleTo< SparsityCsr< ValueType, IndexType > > | |
virtual void | convert_to (result_type *result) const=0 |
Converts the implementer to an object of type result_type. More... | |
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. More... | |
void | move_to (ptr_param< result_type > result) |
Public Member Functions inherited from gko::DiagonalExtractable< ValueType > | |
std::unique_ptr< LinOp > | extract_diagonal_linop () const override |
Extracts the diagonal entries of the matrix into a vector. More... | |
Public Member Functions inherited from gko::ReadableFromMatrixData< ValueType, IndexType > | |
void | read (const matrix_assembly_data< ValueType, IndexType > &data) |
Reads a matrix from a matrix_assembly_data structure. More... | |
Public Member Functions inherited from gko::EnableAbsoluteComputation< remove_complex< Csr< ValueType, IndexType > > > | |
std::unique_ptr< LinOp > | compute_absolute_linop () const override |
Gets the absolute LinOp. More... | |
Public Member Functions inherited from gko::ScaledIdentityAddable | |
void | add_scaled_identity (ptr_param< const LinOp > const a, ptr_param< const LinOp > const b) |
Scales this and adds another scalar times the identity to it. More... | |
Static Public Member Functions | |
static std::unique_ptr< Csr > | create (std::shared_ptr< const Executor > exec, std::shared_ptr< strategy_type > strategy) |
Creates an uninitialized CSR matrix of the specified size. More... | |
static std::unique_ptr< Csr > | create (std::shared_ptr< const Executor > exec, const dim< 2 > &size={}, size_type num_nonzeros={}, std::shared_ptr< strategy_type > strategy=nullptr) |
Creates an uninitialized CSR matrix of the specified size. More... | |
static std::unique_ptr< Csr > | create (std::shared_ptr< const Executor > exec, const dim< 2 > &size, array< value_type > values, array< index_type > col_idxs, array< index_type > row_ptrs, std::shared_ptr< strategy_type > strategy=nullptr) |
Creates a CSR matrix from already allocated (and initialized) row pointer, column index and value arrays. More... | |
template<typename InputValueType , typename InputColumnIndexType , typename InputRowPtrType > | |
static std::unique_ptr< Csr > | create (std::shared_ptr< const Executor > exec, const dim< 2 > &size, std::initializer_list< InputValueType > values, std::initializer_list< InputColumnIndexType > col_idxs, std::initializer_list< InputRowPtrType > row_ptrs) |
create(std::shared_ptr<const Executor>, More... | |
static std::unique_ptr< const Csr > | 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, gko::detail::const_array_view< IndexType > &&row_ptrs, std::shared_ptr< strategy_type > strategy=nullptr) |
Creates a constant (immutable) Csr matrix from a set of constant arrays. More... | |
CSR is a matrix format which stores only the nonzero coefficients by compressing each row of the matrix (compressed sparse row format).
The nonzero elements are stored in a 1D array row-wise, and accompanied with a row pointer array which stores the starting index of each row. An additional column index array is used to identify the column of each nonzero element.
The Csr LinOp supports different operations:
Both the SpGEMM and SpGEAM operation require the input matrices to be sorted by column index, otherwise the algorithms will produce incorrect results.
ValueType | precision of matrix elements |
IndexType | precision of matrix indexes |
gko::matrix::Csr< ValueType, IndexType >::Csr | ( | const Csr< ValueType, IndexType > & | ) |
Copy-constructs a Csr matrix.
Inherits executor, strategy and data.
gko::matrix::Csr< ValueType, IndexType >::Csr | ( | Csr< ValueType, IndexType > && | ) |
|
overridevirtual |
Returns a LinOp representing the column permutation of the Permutable object.
In the resulting LinOp, the column i
contains the input column perm[i]
.
From the linear algebra perspective, with , this represents the operation .
permutation_indices | the array of indices containing the permutation order perm . |
Implements gko::Permutable< IndexType >.
|
overridevirtual |
Gets the AbsoluteLinOp.
Implements gko::EnableAbsoluteComputation< remove_complex< Csr< ValueType, IndexType > > >.
|
overridevirtual |
Returns a LinOp representing the conjugate transpose of the Transposable object.
Implements gko::Transposable.
|
static |
Creates a CSR matrix from already allocated (and initialized) row pointer, 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 |
row_ptrs | array of row pointers |
strategy | the strategy the matrix uses for SpMV operations |
row_ptrs
, col_idxs
or values
is not an rvalue, not an array of IndexType, 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 |
create(std::shared_ptr<const Executor>,
create(std::shared_ptr<const Executor>, const dim<2>&, array<value_type>, array<index_type>, array<index_type>)
References gko::matrix::Csr< ValueType, IndexType >::create().
|
static |
Creates an uninitialized CSR matrix of the specified size.
exec | Executor associated to the matrix |
size | size of the matrix |
num_nonzeros | number of nonzeros |
strategy | the strategy of CSR, or the default strategy if set to nullptr |
|
static |
Creates an uninitialized CSR matrix of the specified size.
exec | Executor associated to the matrix |
strategy | the strategy of CSR |
Referenced by gko::matrix::Csr< ValueType, IndexType >::create().
|
static |
Creates a constant (immutable) Csr 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 |
row_ptrs | the row pointer array of the matrix |
strategy | the strategy the matrix uses for SpMV operations |
std::unique_ptr<Csr<ValueType, IndexType> > gko::matrix::Csr< ValueType, IndexType >::create_submatrix | ( | const index_set< IndexType > & | row_index_set, |
const index_set< IndexType > & | column_index_set | ||
) | const |
Creates a submatrix from this Csr matrix given row and column index_set objects.
row_index_set | the row index set containing the set of rows to be in the submatrix. |
column_index_set | the col index set containing the set of columns to be in the submatrix. |
std::unique_ptr<Csr<ValueType, IndexType> > gko::matrix::Csr< ValueType, IndexType >::create_submatrix | ( | const span & | row_span, |
const span & | column_span | ||
) | const |
Creates a submatrix from this Csr matrix given row and column spans.
row_span | the row span containing the contiguous set of rows to be in the submatrix. |
column_span | the column span containing the contiguous set of columns to be in the submatrix. |
|
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().
|
inlinenoexcept |
Returns the column indexes of the matrix.
References gko::array< ValueType >::get_const_data().
|
inlinenoexcept |
Returns the row pointers of the matrix.
References gko::array< ValueType >::get_const_data().
|
inlinenoexcept |
Returns the starting rows.
References gko::array< ValueType >::get_const_data().
|
inlinenoexcept |
Returns the values of the matrix.
References gko::array< ValueType >::get_const_data().
|
inlinenoexcept |
Returns the number of the srow stored elements (involved warps)
References gko::array< ValueType >::get_size().
|
inlinenoexcept |
Returns the number of elements explicitly stored in the matrix.
References gko::array< ValueType >::get_size().
|
inlinenoexcept |
Returns the row pointers of the matrix.
References gko::array< ValueType >::get_data().
|
inlinenoexcept |
Returns the starting rows.
References gko::array< ValueType >::get_data().
|
inlinenoexcept |
Returns the strategy.
|
inlinenoexcept |
Returns the values of the matrix.
References gko::array< ValueType >::get_data().
|
inline |
Scales the matrix with the inverse of a scalar.
alpha | The entire matrix is scaled by 1 / alpha. alpha has to be a 1x1 Dense matrix. |
References gko::PolymorphicObject::get_executor(), and gko::make_temporary_clone().
|
overridevirtual |
Returns a LinOp representing the row permutation of the inverse permuted object.
In the resulting LinOp, the column perm[i]
contains the input column i
.
From the linear algebra perspective, with , this represents the operation .
permutation_indices | the array of indices containing the permutation order perm . |
Implements gko::Permutable< IndexType >.
|
overridevirtual |
Returns a LinOp representing the symmetric inverse row and column permutation of the Permutable object.
In the resulting LinOp, the entry at location (perm[i],perm[j])
contains the input value (i,j)
.
From the linear algebra perspective, with , this represents the operation .
permutation_indices | the array of indices containing the permutation order. |
Reimplemented from gko::Permutable< IndexType >.
|
overridevirtual |
Returns a LinOp representing the row permutation of the inverse permuted object.
In the resulting LinOp, the row perm[i]
contains the input row i
.
From the linear algebra perspective, with , this represents the operation .
permutation_indices | the array of indices containing the permutation order perm . |
Implements gko::Permutable< IndexType >.
Csr& gko::matrix::Csr< ValueType, IndexType >::operator= | ( | const Csr< ValueType, IndexType > & | ) |
Copy-assigns a Csr matrix.
Preserves executor, copies everything else.
Csr& gko::matrix::Csr< ValueType, IndexType >::operator= | ( | Csr< ValueType, IndexType > && | ) |
|
overridevirtual |
Returns a LinOp representing the symmetric row and column permutation of the Permutable object.
In the resulting LinOp, the entry at location (i,j)
contains the input value (perm[i],perm[j])
.
From the linear algebra perspective, with , this represents the operation .
permutation_indices | the array of indices containing the permutation order. |
Reimplemented from gko::Permutable< IndexType >.
std::unique_ptr<Csr> gko::matrix::Csr< ValueType, IndexType >::permute | ( | ptr_param< const Permutation< index_type >> | permutation, |
permute_mode | mode = permute_mode::symmetric |
||
) | const |
Creates a permuted copy of this matrix with the given permutation .
By default, this computes a symmetric permutation (permute_mode::symmetric). For the effect of the different permutation modes, see permute_mode
permutation | The input permutation. |
mode | The permutation mode. If permute_mode::inverse is set, we use the inverse permutation instead of . If permute_mode::rows is set, the rows will be permuted. If permute_mode::columns is set, the columns will be permuted. |
std::unique_ptr<Csr> gko::matrix::Csr< ValueType, IndexType >::permute | ( | ptr_param< const Permutation< index_type >> | row_permutation, |
ptr_param< const Permutation< index_type >> | column_permutation, | ||
bool | invert = false |
||
) | const |
Creates a non-symmetrically permuted copy of this matrix with the given row and column permutations and .
The operation will compute , or if invert
is false
, and , or if invert
is true
.
row_permutation | The permutation to apply to the rows |
column_permutation | The permutation to apply to the columns |
invert | If set to false , uses the input permutations, otherwise uses their inverses |
|
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 >.
|
overridevirtual |
Returns a LinOp representing the row permutation of the Permutable object.
In the resulting LinOp, the row i
contains the input row perm[i]
.
From the linear algebra perspective, with , this represents the operation .
permutation_indices | the array of indices containing the permutation order. |
Implements gko::Permutable< IndexType >.
|
inline |
Scales the matrix with a scalar.
alpha | The entire matrix is scaled by alpha. alpha has to be a 1x1 Dense matrix. |
References gko::PolymorphicObject::get_executor(), and gko::make_temporary_clone().
std::unique_ptr<Csr> gko::matrix::Csr< ValueType, IndexType >::scale_permute | ( | ptr_param< const ScaledPermutation< value_type, index_type >> | permutation, |
permute_mode | = permute_mode::symmetric |
||
) | const |
Creates a scaled and permuted copy of this matrix.
For an explanation of the permutation modes, see permute(ptr_param<const Permutation<index_type>>, permute_mode)
permutation | The scaled permutation. |
mode | The permutation mode. |
std::unique_ptr<Csr> gko::matrix::Csr< ValueType, IndexType >::scale_permute | ( | ptr_param< const ScaledPermutation< value_type, index_type >> | row_permutation, |
ptr_param< const ScaledPermutation< value_type, index_type >> | column_permutation, | ||
bool | invert = false |
||
) | const |
Creates a scaled and permuted copy of this matrix.
For an explanation of the parameters, see permute(ptr_param<const Permutation<index_type>>, ptr_param<const Permutation<index_type>>, permute_mode)
row_permutation | The scaled row permutation. |
column_permutation | The scaled column permutation. |
invert | If set to false , uses the input permutations, otherwise uses their inverses |
|
inline |
Set the strategy.
strategy | the csr strategy |
|
overridevirtual |
Returns a LinOp representing the transpose of the Transposable object.
Implements gko::Transposable.
|
overridevirtual |
Writes a matrix to a matrix_data structure.
data | the matrix_data structure |
Implements gko::WritableToMatrixData< ValueType, IndexType >.