Ginkgo  Generated from pipelines/2111512875 branch based on develop. Ginkgo version 1.11.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::matrix::Csr< ValueType, IndexType > Class Template Reference

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>

Inheritance diagram for gko::matrix::Csr< ValueType, IndexType >:
[legend]
Collaboration diagram for gko::matrix::Csr< ValueType, IndexType >:
[legend]

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...
 
struct  permuting_reuse_info
 A struct describing a transformation of the matrix that reorders the values of the matrix into the transformed matrix. 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< LinOptranspose () const override
 Returns a LinOp representing the transpose of the Transposable object. More...
 
std::unique_ptr< LinOpconj_transpose () const override
 Returns a LinOp representing the conjugate transpose of the Transposable object. More...
 
std::pair< std::unique_ptr< Csr >, permuting_reuse_infotranspose_reuse () const
 Computes the necessary data to update a transposed matrix from its original matrix. More...
 
std::unique_ptr< Csrpermute (ptr_param< const Permutation< index_type >> permutation, permute_mode mode=permute_mode::symmetric) const
 Creates a permuted copy $A'$ of this matrix $A$ with the given permutation $P$. More...
 
std::unique_ptr< Csrpermute (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 $A'$ of this matrix $A$ with the given row and column permutations $P$ and $Q$. More...
 
std::pair< std::unique_ptr< Csr >, permuting_reuse_infopermute_reuse (ptr_param< const Permutation< index_type >> permutation, permute_mode mode=permute_mode::symmetric) const
 Computes the operations necessary to propagate changed values from a matrix A to a permuted matrix. More...
 
std::pair< std::unique_ptr< Csr >, permuting_reuse_infopermute_reuse (ptr_param< const Permutation< index_type >> row_permutation, ptr_param< const Permutation< index_type >> column_permutation, bool invert=false) const
 Computes the operations necessary to propagate changed values from a matrix A to a permuted matrix. More...
 
std::unique_ptr< Csrscale_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< Csrscale_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< LinOppermute (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< LinOpinverse_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< LinOprow_permute (const array< IndexType > *permutation_indices) const override
 Returns a LinOp representing the row permutation of the Permutable object. More...
 
std::unique_ptr< LinOpcolumn_permute (const array< IndexType > *permutation_indices) const override
 Returns a LinOp representing the column permutation of the Permutable object. More...
 
std::unique_ptr< LinOpinverse_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< LinOpinverse_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...
 
std::unique_ptr< Dense< ValueType > > create_value_view ()
 Creates a Dense view of the value array of this matrix as a column vector of dimensions nnz x 1.
 
std::unique_ptr< const Dense< ValueType > > create_const_value_view () const
 Creates a const Dense view of the value array of this matrix as a column vector of dimensions nnz x 1.
 
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_typeget_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...
 
Csroperator= (const Csr &)
 Copy-assigns a Csr matrix. More...
 
Csroperator= (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
LinOpapply (ptr_param< const LinOp > b, ptr_param< LinOp > x)
 Applies a linear operator to a vector (or a sequence of vectors). More...
 
const LinOpapply (ptr_param< const LinOp > b, ptr_param< LinOp > x) const
 
LinOpapply (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 LinOpapply (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...
 
LinOpoperator= (const LinOp &)=default
 Copy-assigns a LinOp. More...
 
LinOpoperator= (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< LinOpcreate_default (std::shared_ptr< const Executor > exec) const
 
std::unique_ptr< LinOpcreate_default () const
 
std::unique_ptr< LinOpclone (std::shared_ptr< const Executor > exec) const
 
std::unique_ptr< LinOpclone () const
 
LinOpcopy_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)
 
LinOpcopy_from (const std::shared_ptr< const PolymorphicObject > &other)
 
LinOpmove_from (ptr_param< PolymorphicObject > other)
 
LinOpclear ()
 
- Public Member Functions inherited from gko::PolymorphicObject
PolymorphicObjectoperator= (const PolymorphicObject &)
 
std::unique_ptr< PolymorphicObjectcreate_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< PolymorphicObjectcreate_default () const
 Creates a new "default" object of the same dynamic type as this object. More...
 
std::unique_ptr< PolymorphicObjectclone (std::shared_ptr< const Executor > exec) const
 Creates a clone of the object. More...
 
std::unique_ptr< PolymorphicObjectclone () const
 Creates a clone of the object. More...
 
PolymorphicObjectcopy_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...
 
PolymorphicObjectcopy_from (const std::shared_ptr< const PolymorphicObject > &other)
 Copies another object into this object. More...
 
PolymorphicObjectmove_from (ptr_param< PolymorphicObject > other)
 Moves another object into this object. More...
 
PolymorphicObjectclear ()
 Transforms the object into its default state. More...
 
std::shared_ptr< const Executorget_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< LinOpextract_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< LinOpcompute_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< Csrcreate (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< Csrcreate (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< Csrcreate (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< Csrcreate (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 Csrcreate_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...
 

Friends

class EnablePolymorphicObject< Csr, LinOp >
 
class Coo< ValueType, IndexType >
 
class Dense< ValueType >
 
class Diagonal< ValueType >
 
class Ell< ValueType, IndexType >
 
class Hybrid< ValueType, IndexType >
 
class Sellp< ValueType, IndexType >
 
class SparsityCsr< ValueType, IndexType >
 
class Fbcsr< ValueType, IndexType >
 
class CsrBuilder< ValueType, IndexType >
 
class Csr< to_complex< ValueType >, IndexType >
 
class Csr< previous_precision< ValueType >, IndexType >
 

Detailed Description

template<typename ValueType = default_precision, typename IndexType = int32>
class gko::matrix::Csr< ValueType, IndexType >

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:

matrix::Csr *A, *B, *C; // matrices
matrix::Dense *b, *x; // vectors tall-and-skinny matrices
matrix::Dense *alpha, *beta; // scalars of dimension 1x1
matrix::Identity *I; // identity matrix
// Applying to Dense matrices computes an SpMV/SpMM product
A->apply(b, x) // x = A*b
A->apply(alpha, b, beta, x) // x = alpha*A*b + beta*x
// Applying to Csr matrices computes a SpGEMM product of two sparse matrices
A->apply(B, C) // C = A*B
A->apply(alpha, B, beta, C) // C = alpha*A*B + beta*C
// Applying to an Identity matrix computes a SpGEAM sparse matrix addition
A->apply(alpha, I, beta, B) // B = alpha*A + beta*B

Both the SpGEMM and SpGEAM operation require the input matrices to be sorted by column index, otherwise the algorithms will produce incorrect results.

Template Parameters
ValueTypeprecision of matrix elements
IndexTypeprecision of matrix indexes

Constructor & Destructor Documentation

◆ Csr() [1/2]

template<typename ValueType = default_precision, typename IndexType = int32>
gko::matrix::Csr< ValueType, IndexType >::Csr ( const Csr< ValueType, IndexType > &  )

Copy-constructs a Csr matrix.

Inherits executor, strategy and data.

◆ Csr() [2/2]

template<typename ValueType = default_precision, typename IndexType = int32>
gko::matrix::Csr< ValueType, IndexType >::Csr ( Csr< ValueType, IndexType > &&  )

Move-constructs a Csr matrix.

Inherits executor and strategy, moves the data and leaves the moved-from object in an empty state (0x0 LinOp with unchanged executor and strategy, no nonzeros and valid row pointers).

Member Function Documentation

◆ column_permute()

template<typename ValueType = default_precision, typename IndexType = int32>
std::unique_ptr<LinOp> gko::matrix::Csr< ValueType, IndexType >::column_permute ( const array< IndexType > *  permutation_indices) const
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 $P_{ij} = \delta_{i \pi(i)}$, this represents the operation $A P^T$.

Parameters
permutation_indicesthe array of indices containing the permutation order perm.
Returns
a pointer to the new column permuted object

Implements gko::Permutable< IndexType >.

◆ compute_absolute()

template<typename ValueType = default_precision, typename IndexType = int32>
std::unique_ptr<absolute_type> gko::matrix::Csr< ValueType, IndexType >::compute_absolute ( ) const
overridevirtual

Gets the AbsoluteLinOp.

Returns
a pointer to the new absolute object

Implements gko::EnableAbsoluteComputation< remove_complex< Csr< ValueType, IndexType > > >.

◆ conj_transpose()

template<typename ValueType = default_precision, typename IndexType = int32>
std::unique_ptr<LinOp> gko::matrix::Csr< ValueType, IndexType >::conj_transpose ( ) const
overridevirtual

Returns a LinOp representing the conjugate transpose of the Transposable object.

Returns
a pointer to the new conjugate transposed object

Implements gko::Transposable.

◆ create() [1/4]

template<typename ValueType = default_precision, typename IndexType = int32>
static std::unique_ptr<Csr> gko::matrix::Csr< ValueType, IndexType >::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 
)
static

Creates a CSR matrix from already allocated (and initialized) row pointer, column index and value arrays.

Parameters
execExecutor associated to the matrix
sizesize of the matrix
valuesarray of matrix values
col_idxsarray of column indexes
row_ptrsarray of row pointers
strategythe strategy the matrix uses for SpMV operations
Note
If one of 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.
Returns
A smart pointer to the newly created matrix.

◆ create() [2/4]

template<typename ValueType = default_precision, typename IndexType = int32>
template<typename InputValueType , typename InputColumnIndexType , typename InputRowPtrType >
static std::unique_ptr<Csr> gko::matrix::Csr< ValueType, IndexType >::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 
)
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().

◆ create() [3/4]

template<typename ValueType = default_precision, typename IndexType = int32>
static std::unique_ptr<Csr> gko::matrix::Csr< ValueType, IndexType >::create ( std::shared_ptr< const Executor exec,
const dim< 2 > &  size = {},
size_type  num_nonzeros = {},
std::shared_ptr< strategy_type strategy = nullptr 
)
static

Creates an uninitialized CSR matrix of the specified size.

Parameters
execExecutor associated to the matrix
sizesize of the matrix
num_nonzerosnumber of nonzeros
strategythe strategy of CSR, or the default strategy if set to nullptr
Returns
A smart pointer to the newly created matrix.

◆ create() [4/4]

template<typename ValueType = default_precision, typename IndexType = int32>
static std::unique_ptr<Csr> gko::matrix::Csr< ValueType, IndexType >::create ( std::shared_ptr< const Executor exec,
std::shared_ptr< strategy_type strategy 
)
static

Creates an uninitialized CSR matrix of the specified size.

Parameters
execExecutor associated to the matrix
strategythe strategy of CSR
Returns
A smart pointer to the newly created matrix.

Referenced by gko::matrix::Csr< ValueType, IndexType >::create().

◆ create_const()

template<typename ValueType = default_precision, typename IndexType = int32>
static std::unique_ptr<const Csr> gko::matrix::Csr< ValueType, IndexType >::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 
)
static

Creates a constant (immutable) Csr matrix from a set of constant arrays.

Parameters
execthe executor to create the matrix on
sizethe dimensions of the matrix
valuesthe value array of the matrix
col_idxsthe column index array of the matrix
row_ptrsthe row pointer array of the matrix
strategythe strategy the matrix uses for SpMV operations
Returns
A smart pointer to the constant matrix wrapping the input arrays (if they reside on the same executor as the matrix) or a copy of these arrays on the correct executor.
A smart pointer to the newly created matrix.

◆ create_submatrix() [1/2]

template<typename ValueType = default_precision, typename IndexType = int32>
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.

Parameters
row_index_setthe row index set containing the set of rows to be in the submatrix.
column_index_setthe col index set containing the set of columns to be in the submatrix.
Returns
A new CSR matrix with the elements that belong to the row and columns of this matrix as specified by the index sets.
Note
This is not a view but creates a new, separate CSR matrix.

◆ create_submatrix() [2/2]

template<typename ValueType = default_precision, typename IndexType = int32>
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.

Parameters
row_spanthe row span containing the contiguous set of rows to be in the submatrix.
column_spanthe column span containing the contiguous set of columns to be in the submatrix.
Returns
A new CSR matrix with the elements that belong to the row and columns of this matrix as specified by the index sets.
Note
This is not a view but creates a new, separate CSR matrix.

◆ extract_diagonal()

template<typename ValueType = default_precision, typename IndexType = int32>
std::unique_ptr<Diagonal<ValueType> > gko::matrix::Csr< ValueType, IndexType >::extract_diagonal ( ) const
overridevirtual

Extracts the diagonal entries of the matrix into a vector.

Parameters
diagthe vector into which the diagonal will be written

Implements gko::DiagonalExtractable< ValueType >.

◆ get_col_idxs()

template<typename ValueType = default_precision, typename IndexType = int32>
index_type* gko::matrix::Csr< ValueType, IndexType >::get_col_idxs ( )
inlinenoexcept

Returns the column indexes of the matrix.

Returns
the column indexes of the matrix.

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

◆ get_const_col_idxs()

template<typename ValueType = default_precision, typename IndexType = int32>
const index_type* gko::matrix::Csr< ValueType, IndexType >::get_const_col_idxs ( ) const
inlinenoexcept

Returns the column indexes of the matrix.

Returns
the column indexes of the matrix.
Note
This is the constant version of the function, which can be significantly more memory efficient than the non-constant version, so always prefer this version.

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

◆ get_const_row_ptrs()

template<typename ValueType = default_precision, typename IndexType = int32>
const index_type* gko::matrix::Csr< ValueType, IndexType >::get_const_row_ptrs ( ) const
inlinenoexcept

Returns the row pointers of the matrix.

Returns
the row pointers of the matrix.
Note
This is the constant version of the function, which can be significantly more memory efficient than the non-constant version, so always prefer this version.

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

◆ get_const_srow()

template<typename ValueType = default_precision, typename IndexType = int32>
const index_type* gko::matrix::Csr< ValueType, IndexType >::get_const_srow ( ) const
inlinenoexcept

Returns the starting rows.

Returns
the starting rows.
Note
This is the constant version of the function, which can be significantly more memory efficient than the non-constant version, so always prefer this version.

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

◆ get_const_values()

template<typename ValueType = default_precision, typename IndexType = int32>
const value_type* gko::matrix::Csr< ValueType, IndexType >::get_const_values ( ) const
inlinenoexcept

Returns the values of the matrix.

Returns
the values of the matrix.
Note
This is the constant version of the function, which can be significantly more memory efficient than the non-constant version, so always prefer this version.

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

◆ get_num_srow_elements()

template<typename ValueType = default_precision, typename IndexType = int32>
size_type gko::matrix::Csr< ValueType, IndexType >::get_num_srow_elements ( ) const
inlinenoexcept

Returns the number of the srow stored elements (involved warps)

Returns
the number of the srow stored elements (involved warps)

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

◆ get_num_stored_elements()

template<typename ValueType = default_precision, typename IndexType = int32>
size_type gko::matrix::Csr< ValueType, IndexType >::get_num_stored_elements ( ) const
inlinenoexcept

Returns the number of elements explicitly stored in the matrix.

Returns
the number of elements explicitly stored in the matrix

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

◆ get_row_ptrs()

template<typename ValueType = default_precision, typename IndexType = int32>
index_type* gko::matrix::Csr< ValueType, IndexType >::get_row_ptrs ( )
inlinenoexcept

Returns the row pointers of the matrix.

Returns
the row pointers of the matrix.

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

◆ get_srow()

template<typename ValueType = default_precision, typename IndexType = int32>
index_type* gko::matrix::Csr< ValueType, IndexType >::get_srow ( )
inlinenoexcept

Returns the starting rows.

Returns
the starting rows.

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

◆ get_strategy()

template<typename ValueType = default_precision, typename IndexType = int32>
std::shared_ptr<strategy_type> gko::matrix::Csr< ValueType, IndexType >::get_strategy ( ) const
inlinenoexcept

Returns the strategy.

Returns
the strategy

◆ get_values()

template<typename ValueType = default_precision, typename IndexType = int32>
value_type* gko::matrix::Csr< ValueType, IndexType >::get_values ( )
inlinenoexcept

Returns the values of the matrix.

Returns
the values of the matrix.

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

◆ inv_scale()

template<typename ValueType = default_precision, typename IndexType = int32>
void gko::matrix::Csr< ValueType, IndexType >::inv_scale ( ptr_param< const LinOp alpha)
inline

Scales the matrix with the inverse of a scalar.

Parameters
alphaThe 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().

◆ inverse_column_permute()

template<typename ValueType = default_precision, typename IndexType = int32>
std::unique_ptr<LinOp> gko::matrix::Csr< ValueType, IndexType >::inverse_column_permute ( const array< IndexType > *  permutation_indices) const
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 $P_{ij} = \delta_{i \pi(i)}$, this represents the operation $A P^{-T}$.

Parameters
permutation_indicesthe array of indices containing the permutation order perm.
Returns
a pointer to the new inverse permuted object

Implements gko::Permutable< IndexType >.

◆ inverse_permute()

template<typename ValueType = default_precision, typename IndexType = int32>
std::unique_ptr<LinOp> gko::matrix::Csr< ValueType, IndexType >::inverse_permute ( const array< IndexType > *  permutation_indices) const
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 $P_{ij} = \delta_{i \pi(i)}$, this represents the operation $P^{-1} A P^{-T}$.

Parameters
permutation_indicesthe array of indices containing the permutation order.
Returns
a pointer to the new permuted object

Reimplemented from gko::Permutable< IndexType >.

◆ inverse_row_permute()

template<typename ValueType = default_precision, typename IndexType = int32>
std::unique_ptr<LinOp> gko::matrix::Csr< ValueType, IndexType >::inverse_row_permute ( const array< IndexType > *  permutation_indices) const
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 $P_{ij} = \delta_{i \pi(i)}$, this represents the operation $P^{-1} A$.

Parameters
permutation_indicesthe array of indices containing the permutation order perm.
Returns
a pointer to the new inverse permuted object

Implements gko::Permutable< IndexType >.

◆ operator=() [1/2]

template<typename ValueType = default_precision, typename IndexType = int32>
Csr& gko::matrix::Csr< ValueType, IndexType >::operator= ( const Csr< ValueType, IndexType > &  )

Copy-assigns a Csr matrix.

Preserves executor, copies everything else.

◆ operator=() [2/2]

template<typename ValueType = default_precision, typename IndexType = int32>
Csr& gko::matrix::Csr< ValueType, IndexType >::operator= ( Csr< ValueType, IndexType > &&  )

Move-assigns a Csr matrix.

Preserves executor, moves the data and leaves the moved-from object in an empty state (0x0 LinOp with unchanged executor and strategy, no nonzeros and valid row pointers).

◆ permute() [1/3]

template<typename ValueType = default_precision, typename IndexType = int32>
std::unique_ptr<LinOp> gko::matrix::Csr< ValueType, IndexType >::permute ( const array< IndexType > *  permutation_indices) const
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 $P_{ij} = \delta_{i \pi(i)}$, this represents the operation $P A P^T$.

Parameters
permutation_indicesthe array of indices containing the permutation order.
Returns
a pointer to the new permuted object

Reimplemented from gko::Permutable< IndexType >.

◆ permute() [2/3]

template<typename ValueType = default_precision, typename IndexType = int32>
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 $A'$ of this matrix $A$ with the given permutation $P$.

By default, this computes a symmetric permutation (permute_mode::symmetric). For the effect of the different permutation modes, see permute_mode

Parameters
permutationThe input permutation.
modeThe permutation mode. If permute_mode::inverse is set, we use the inverse permutation $P^{-1}$ instead of $P$. If permute_mode::rows is set, the rows will be permuted. If permute_mode::columns is set, the columns will be permuted.
Returns
The permuted matrix.

◆ permute() [3/3]

template<typename ValueType = default_precision, typename IndexType = int32>
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 $A'$ of this matrix $A$ with the given row and column permutations $P$ and $Q$.

The operation will compute $A'(i, j) = A(p[i], q[j])$, or $A' = P A Q^T$ if invert is false, and $A'(p[i], q[j]) = A(i,j)$, or $A' = P^{-1} A Q^{-T}$ if invert is true.

Parameters
row_permutationThe permutation $P$ to apply to the rows
column_permutationThe permutation $Q$ to apply to the columns
invertIf set to false, uses the input permutations, otherwise uses their inverses $P^{-1}, Q^{-1}$
Returns
The permuted matrix.

◆ permute_reuse() [1/2]

template<typename ValueType = default_precision, typename IndexType = int32>
std::pair<std::unique_ptr<Csr>, permuting_reuse_info> gko::matrix::Csr< ValueType, IndexType >::permute_reuse ( ptr_param< const Permutation< index_type >>  permutation,
permute_mode  mode = permute_mode::symmetric 
) const

Computes the operations necessary to propagate changed values from a matrix A to a permuted matrix.

The semantics of this function match those of permute(ptr_param<const Permutation<index_type>>, permute_mode). Updating values works as follows:

auto [permuted, reuse] = matrix->permute_reuse(permutation, mode);
change_values(matrix);
reuse->update_values(matrix, permuted);
Parameters
permutationThe input permutation.
modeThe permutation mode. If permute_mode::inverse is set, we use the inverse permutation $P^{-1}$ instead of $P$. If permute_mode::rows is set, the rows will be permuted. If permute_mode::columns is set, the columns will be permuted.
Returns
an std::pair consisting of the permuted matrix and the reuse info that can be used to update values in the permuted matrix.

◆ permute_reuse() [2/2]

template<typename ValueType = default_precision, typename IndexType = int32>
std::pair<std::unique_ptr<Csr>, permuting_reuse_info> gko::matrix::Csr< ValueType, IndexType >::permute_reuse ( ptr_param< const Permutation< index_type >>  row_permutation,
ptr_param< const Permutation< index_type >>  column_permutation,
bool  invert = false 
) const

Computes the operations necessary to propagate changed values from a matrix A to a permuted matrix.

The semantics of this function match those of permute(ptr_param<const Permutation<index_type>>, ptr_param<const Permutation<index_type>>, bool). Updating values works as follows:

auto [permuted, reuse] = matrix->permute_reuse(row_perm, col_perm, inv);
change_values(matrix);
reuse->update_values(matrix, permuted);
Parameters
row_permutationThe permutation $P$ to apply to the rows
column_permutationThe permutation $Q$ to apply to the columns
invertIf set to false, uses the input permutations, otherwise uses their inverses $P^{-1}, Q^{-1}$
Returns
an std::pair consisting of the permuted matrix and the reuse info that can be used to update values in the permuted matrix.

◆ read() [1/3]

template<typename ValueType = default_precision, typename IndexType = int32>
void gko::matrix::Csr< ValueType, IndexType >::read ( const device_mat_data data)
overridevirtual

Reads a matrix from a device_matrix_data structure.

Parameters
datathe device_matrix_data structure.

Reimplemented from gko::ReadableFromMatrixData< ValueType, IndexType >.

◆ read() [2/3]

template<typename ValueType = default_precision, typename IndexType = int32>
void gko::matrix::Csr< ValueType, IndexType >::read ( const mat_data data)
overridevirtual

Reads a matrix from a matrix_data structure.

Parameters
datathe matrix_data structure

Implements gko::ReadableFromMatrixData< ValueType, IndexType >.

◆ read() [3/3]

template<typename ValueType = default_precision, typename IndexType = int32>
void gko::matrix::Csr< ValueType, IndexType >::read ( device_mat_data &&  data)
overridevirtual

Reads a matrix from a device_matrix_data structure.

The structure may be emptied by this function.

Parameters
datathe device_matrix_data structure.

Reimplemented from gko::ReadableFromMatrixData< ValueType, IndexType >.

◆ row_permute()

template<typename ValueType = default_precision, typename IndexType = int32>
std::unique_ptr<LinOp> gko::matrix::Csr< ValueType, IndexType >::row_permute ( const array< IndexType > *  permutation_indices) const
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 $P_{ij} = \delta_{i \pi(i)}$, this represents the operation $P A$.

Parameters
permutation_indicesthe array of indices containing the permutation order.
Returns
a pointer to the new permuted object

Implements gko::Permutable< IndexType >.

◆ scale()

template<typename ValueType = default_precision, typename IndexType = int32>
void gko::matrix::Csr< ValueType, IndexType >::scale ( ptr_param< const LinOp alpha)
inline

Scales the matrix with a scalar.

Parameters
alphaThe entire matrix is scaled by alpha. alpha has to be a 1x1 Dense matrix.

References gko::PolymorphicObject::get_executor(), and gko::make_temporary_clone().

◆ scale_permute() [1/2]

template<typename ValueType = default_precision, typename IndexType = int32>
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)

Parameters
permutationThe scaled permutation.
modeThe permutation mode.
Returns
The permuted matrix.

◆ scale_permute() [2/2]

template<typename ValueType = default_precision, typename IndexType = int32>
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)

Parameters
row_permutationThe scaled row permutation.
column_permutationThe scaled column permutation.
invertIf set to false, uses the input permutations, otherwise uses their inverses $P^{-1}, Q^{-1}$
Returns
The permuted matrix.

◆ set_strategy()

template<typename ValueType = default_precision, typename IndexType = int32>
void gko::matrix::Csr< ValueType, IndexType >::set_strategy ( std::shared_ptr< strategy_type strategy)
inline

Set the strategy.

Parameters
strategythe csr strategy

◆ transpose()

template<typename ValueType = default_precision, typename IndexType = int32>
std::unique_ptr<LinOp> gko::matrix::Csr< ValueType, IndexType >::transpose ( ) const
overridevirtual

Returns a LinOp representing the transpose of the Transposable object.

Returns
a pointer to the new transposed object

Implements gko::Transposable.

◆ transpose_reuse()

template<typename ValueType = default_precision, typename IndexType = int32>
std::pair<std::unique_ptr<Csr>, permuting_reuse_info> gko::matrix::Csr< ValueType, IndexType >::transpose_reuse ( ) const

Computes the necessary data to update a transposed matrix from its original matrix.

auto [transposed, reuse] = matrix->transpose_reuse();
change_values(matrix);
reuse->update_values(matrix, transposed);
Returns
the reuse info struct that can be used to update values in the transposed matrix.

◆ write()

template<typename ValueType = default_precision, typename IndexType = int32>
void gko::matrix::Csr< ValueType, IndexType >::write ( mat_data data) const
overridevirtual

Writes a matrix to a matrix_data structure.

Parameters
datathe matrix_data structure

Implements gko::WritableToMatrixData< ValueType, IndexType >.


The documentation for this class was generated from the following files:
gko::stop::mode
mode
The mode for the residual norm criterion.
Definition: residual_norm.hpp:38