Ginkgo  Generated from pipelines/1068515030 branch based on master. Ginkgo version 1.7.0
A numerical linear algebra library targeting many-core architectures
Public Types | Public Member Functions | Static Public Member Functions | Friends | List of all members
gko::matrix::Dense< ValueType > Class Template Reference

Dense is a matrix format which explicitly stores all values of the matrix. More...

#include <ginkgo/core/matrix/dense.hpp>

Inheritance diagram for gko::matrix::Dense< ValueType >:
[legend]
Collaboration diagram for gko::matrix::Dense< ValueType >:
[legend]

Public Types

using value_type = ValueType
 
using index_type = int64
 
using transposed_type = Dense< ValueType >
 
using mat_data = matrix_data< ValueType, int64 >
 
using mat_data32 = matrix_data< ValueType, int32 >
 
using device_mat_data = device_matrix_data< ValueType, int64 >
 
using device_mat_data32 = device_matrix_data< ValueType, int32 >
 
using absolute_type = remove_complex< Dense >
 
using real_type = absolute_type
 
using complex_type = to_complex< Dense >
 
using row_major_range = gko::range< gko::accessor::row_major< ValueType, 2 > >
 
- Public Types inherited from gko::EnablePolymorphicAssignment< Dense< ValueType > >
using result_type = Dense< ValueType >
 
- Public Types inherited from gko::ConvertibleTo< Dense< ValueType > >
using result_type = Dense< ValueType >
 
- Public Types inherited from gko::ConvertibleTo< Dense< next_precision< ValueType > > >
using result_type = Dense< next_precision< ValueType > >
 
- Public Types inherited from gko::ConvertibleTo< Coo< ValueType, int32 > >
using result_type = Coo< ValueType, int32 >
 
- Public Types inherited from gko::ConvertibleTo< Coo< ValueType, int64 > >
using result_type = Coo< ValueType, int64 >
 
- Public Types inherited from gko::ConvertibleTo< Csr< ValueType, int32 > >
using result_type = Csr< ValueType, int32 >
 
- Public Types inherited from gko::ConvertibleTo< Csr< ValueType, int64 > >
using result_type = Csr< ValueType, int64 >
 
- Public Types inherited from gko::ConvertibleTo< Ell< ValueType, int32 > >
using result_type = Ell< ValueType, int32 >
 
- Public Types inherited from gko::ConvertibleTo< Ell< ValueType, int64 > >
using result_type = Ell< ValueType, int64 >
 
- Public Types inherited from gko::ConvertibleTo< Fbcsr< ValueType, int32 > >
using result_type = Fbcsr< ValueType, int32 >
 
- Public Types inherited from gko::ConvertibleTo< Fbcsr< ValueType, int64 > >
using result_type = Fbcsr< ValueType, int64 >
 
- Public Types inherited from gko::ConvertibleTo< Hybrid< ValueType, int32 > >
using result_type = Hybrid< ValueType, int32 >
 
- Public Types inherited from gko::ConvertibleTo< Hybrid< ValueType, int64 > >
using result_type = Hybrid< ValueType, int64 >
 
- Public Types inherited from gko::ConvertibleTo< Sellp< ValueType, int32 > >
using result_type = Sellp< ValueType, int32 >
 
- Public Types inherited from gko::ConvertibleTo< Sellp< ValueType, int64 > >
using result_type = Sellp< ValueType, int64 >
 
- Public Types inherited from gko::ConvertibleTo< SparsityCsr< ValueType, int32 > >
using result_type = SparsityCsr< ValueType, int32 >
 
- Public Types inherited from gko::ConvertibleTo< SparsityCsr< ValueType, int64 > >
using result_type = SparsityCsr< ValueType, int64 >
 
- Public Types inherited from gko::DiagonalExtractable< ValueType >
using value_type = ValueType
 
- Public Types inherited from gko::ReadableFromMatrixData< ValueType, int32 >
using value_type = ValueType
 
using index_type = int32
 
- Public Types inherited from gko::ReadableFromMatrixData< ValueType, int64 >
using value_type = ValueType
 
using index_type = int64
 
- Public Types inherited from gko::WritableToMatrixData< ValueType, int32 >
using value_type = ValueType
 
using index_type = int32
 
- Public Types inherited from gko::WritableToMatrixData< ValueType, int64 >
using value_type = ValueType
 
using index_type = int64
 
- Public Types inherited from gko::EnableAbsoluteComputation< remove_complex< Dense< ValueType > > >
using absolute_type = remove_complex< Dense< ValueType > >
 

Public Member Functions

void convert_to (Dense< next_precision< ValueType >> *result) const override
 
void move_to (Dense< next_precision< ValueType >> *result) override
 
void convert_to (Coo< ValueType, int32 > *result) const override
 
void move_to (Coo< ValueType, int32 > *result) override
 
void convert_to (Coo< ValueType, int64 > *result) const override
 
void move_to (Coo< ValueType, int64 > *result) override
 
void convert_to (Csr< ValueType, int32 > *result) const override
 
void move_to (Csr< ValueType, int32 > *result) override
 
void convert_to (Csr< ValueType, int64 > *result) const override
 
void move_to (Csr< ValueType, int64 > *result) override
 
void convert_to (Ell< ValueType, int32 > *result) const override
 
void move_to (Ell< ValueType, int32 > *result) override
 
void convert_to (Ell< ValueType, int64 > *result) const override
 
void move_to (Ell< ValueType, int64 > *result) override
 
void convert_to (Fbcsr< ValueType, int32 > *result) const override
 
void move_to (Fbcsr< ValueType, int32 > *result) override
 
void convert_to (Fbcsr< ValueType, int64 > *result) const override
 
void move_to (Fbcsr< ValueType, int64 > *result) override
 
void convert_to (Hybrid< ValueType, int32 > *result) const override
 
void move_to (Hybrid< ValueType, int32 > *result) override
 
void convert_to (Hybrid< ValueType, int64 > *result) const override
 
void move_to (Hybrid< ValueType, int64 > *result) override
 
void convert_to (Sellp< ValueType, int32 > *result) const override
 
void move_to (Sellp< ValueType, int32 > *result) override
 
void convert_to (Sellp< ValueType, int64 > *result) const override
 
void move_to (Sellp< ValueType, int64 > *result) override
 
void convert_to (SparsityCsr< ValueType, int32 > *result) const override
 
void move_to (SparsityCsr< ValueType, int32 > *result) override
 
void convert_to (SparsityCsr< ValueType, int64 > *result) const override
 
void move_to (SparsityCsr< ValueType, int64 > *result) override
 
void read (const mat_data &data) override
 
void read (const mat_data32 &data) override
 
void read (const device_mat_data &data) override
 
void read (const device_mat_data32 &data) override
 
void read (device_mat_data &&data) override
 
void read (device_mat_data32 &&data) override
 
void write (mat_data &data) const override
 
void write (mat_data32 &data) const override
 
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...
 
void transpose (ptr_param< Dense > output) const
 Writes the transposed matrix into the given output matrix. More...
 
void conj_transpose (ptr_param< Dense > output) const
 Writes the conjugate-transposed matrix into the given output matrix. More...
 
void fill (const ValueType value)
 Fill the dense matrix with a given value. More...
 
std::unique_ptr< Densepermute (ptr_param< const Permutation< int32 >> 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< Densepermute (ptr_param< const Permutation< int64 >> permutation, permute_mode mode=permute_mode::symmetric) const
 
void permute (ptr_param< const Permutation< int32 >> permutation, ptr_param< Dense > output, permute_mode mode) const
 Overload of permute(ptr_param<const Permutation<int32>>, permute_mode) that writes the permuted copy into an existing Dense matrix. More...
 
void permute (ptr_param< const Permutation< int64 >> permutation, ptr_param< Dense > output, permute_mode mode) const
 
std::unique_ptr< Densepermute (ptr_param< const Permutation< int32 >> row_permutation, ptr_param< const Permutation< int32 >> 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::unique_ptr< Densepermute (ptr_param< const Permutation< int64 >> row_permutation, ptr_param< const Permutation< int64 >> column_permutation, bool invert=false) const
 
void permute (ptr_param< const Permutation< int32 >> row_permutation, ptr_param< const Permutation< int32 >> column_permutation, ptr_param< Dense > output, bool invert=false) const
 Overload of permute(ptr_param<const Permutation<int32>>, ptr_param<const Permutation<int32>>, permute_mode) that writes the permuted copy into an existing Dense matrix. More...
 
void permute (ptr_param< const Permutation< int64 >> row_permutation, ptr_param< const Permutation< int64 >> column_permutation, ptr_param< Dense > output, bool invert=false) const
 
std::unique_ptr< Densescale_permute (ptr_param< const ScaledPermutation< value_type, int32 >> permutation, permute_mode mode=permute_mode::symmetric) const
 Creates a scaled and permuted copy of this matrix. More...
 
std::unique_ptr< Densescale_permute (ptr_param< const ScaledPermutation< value_type, int64 >> permutation, permute_mode mode=permute_mode::symmetric) const
 
void scale_permute (ptr_param< const ScaledPermutation< value_type, int32 >> permutation, ptr_param< Dense > output, permute_mode mode) const
 Overload of scale_permute(ptr_param<const ScaledPermutation<value_type, int32>>, permute_mode) that writes the permuted copy into an existing Dense matrix. More...
 
void scale_permute (ptr_param< const ScaledPermutation< value_type, int64 >> permutation, ptr_param< Dense > output, permute_mode mode) const
 
std::unique_ptr< Densescale_permute (ptr_param< const ScaledPermutation< value_type, int32 >> row_permutation, ptr_param< const ScaledPermutation< value_type, int32 >> column_permutation, bool invert=false) const
 Creates a scaled and permuted copy of this matrix. More...
 
std::unique_ptr< Densescale_permute (ptr_param< const ScaledPermutation< value_type, int64 >> row_permutation, ptr_param< const ScaledPermutation< value_type, int64 >> column_permutation, bool invert=false) const
 
void scale_permute (ptr_param< const ScaledPermutation< value_type, int32 >> row_permutation, ptr_param< const ScaledPermutation< value_type, int32 >> column_permutation, ptr_param< Dense > output, bool invert=false) const
 Overload of scale_permute(ptr_param<const ScaledPermutation<value_type, int32>>, ptr_param<const ScaledPermutation<value_type, int32>>, bool) that writes the permuted copy into an existing Dense matrix. More...
 
void scale_permute (ptr_param< const ScaledPermutation< value_type, int64 >> row_permutation, ptr_param< const ScaledPermutation< value_type, int64 >> column_permutation, ptr_param< Dense > output, bool invert=false) const
 
std::unique_ptr< LinOppermute (const array< int32 > *permutation_indices) const override
 Returns a LinOp representing the symmetric row and column permutation of the Permutable object. More...
 
std::unique_ptr< LinOppermute (const array< int64 > *permutation_indices) const override
 Returns a LinOp representing the symmetric row and column permutation of the Permutable object. More...
 
void permute (const array< int32 > *permutation_indices, ptr_param< Dense > output) const
 Writes the symmetrically permuted matrix into the given output matrix. More...
 
void permute (const array< int64 > *permutation_indices, ptr_param< Dense > output) const
 
std::unique_ptr< LinOpinverse_permute (const array< int32 > *permutation_indices) const override
 Returns a LinOp representing the symmetric inverse row and column permutation of the Permutable object. More...
 
std::unique_ptr< LinOpinverse_permute (const array< int64 > *permutation_indices) const override
 Returns a LinOp representing the symmetric inverse row and column permutation of the Permutable object. More...
 
void inverse_permute (const array< int32 > *permutation_indices, ptr_param< Dense > output) const
 Writes the inverse symmetrically permuted matrix into the given output matrix. More...
 
void inverse_permute (const array< int64 > *permutation_indices, ptr_param< Dense > output) const
 
std::unique_ptr< LinOprow_permute (const array< int32 > *permutation_indices) const override
 Returns a LinOp representing the row permutation of the Permutable object. More...
 
std::unique_ptr< LinOprow_permute (const array< int64 > *permutation_indices) const override
 Returns a LinOp representing the row permutation of the Permutable object. More...
 
void row_permute (const array< int32 > *permutation_indices, ptr_param< Dense > output) const
 Writes the row-permuted matrix into the given output matrix. More...
 
void row_permute (const array< int64 > *permutation_indices, ptr_param< Dense > output) const
 
std::unique_ptr< Denserow_gather (const array< int32 > *gather_indices) const
 Create a Dense matrix consisting of the given rows from this matrix. More...
 
std::unique_ptr< Denserow_gather (const array< int64 > *gather_indices) const
 Create a Dense matrix consisting of the given rows from this matrix. More...
 
void row_gather (const array< int32 > *gather_indices, ptr_param< LinOp > row_collection) const
 Copies the given rows from this matrix into row_collection More...
 
void row_gather (const array< int64 > *gather_indices, ptr_param< LinOp > row_collection) const
 
void row_gather (ptr_param< const LinOp > alpha, const array< int32 > *gather_indices, ptr_param< const LinOp > beta, ptr_param< LinOp > row_collection) const
 Copies the given rows from this matrix into row_collection with scaling. More...
 
void row_gather (ptr_param< const LinOp > alpha, const array< int64 > *gather_indices, ptr_param< const LinOp > beta, ptr_param< LinOp > row_collection) const
 
std::unique_ptr< LinOpcolumn_permute (const array< int32 > *permutation_indices) const override
 Returns a LinOp representing the column permutation of the Permutable object. More...
 
std::unique_ptr< LinOpcolumn_permute (const array< int64 > *permutation_indices) const override
 Returns a LinOp representing the column permutation of the Permutable object. More...
 
void column_permute (const array< int32 > *permutation_indices, ptr_param< Dense > output) const
 Writes the column-permuted matrix into the given output matrix. More...
 
void column_permute (const array< int64 > *permutation_indices, ptr_param< Dense > output) const
 
std::unique_ptr< LinOpinverse_row_permute (const array< int32 > *permutation_indices) const override
 Returns a LinOp representing the row permutation of the inverse permuted object. More...
 
std::unique_ptr< LinOpinverse_row_permute (const array< int64 > *permutation_indices) const override
 Returns a LinOp representing the row permutation of the inverse permuted object. More...
 
void inverse_row_permute (const array< int32 > *permutation_indices, ptr_param< Dense > output) const
 Writes the inverse row-permuted matrix into the given output matrix. More...
 
void inverse_row_permute (const array< int64 > *permutation_indices, ptr_param< Dense > output) const
 
std::unique_ptr< LinOpinverse_column_permute (const array< int32 > *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< int64 > *permutation_indices) const override
 Returns a LinOp representing the row permutation of the inverse permuted object. More...
 
void inverse_column_permute (const array< int32 > *permutation_indices, ptr_param< Dense > output) const
 Writes the inverse column-permuted matrix into the given output matrix. More...
 
void inverse_column_permute (const array< int64 > *permutation_indices, ptr_param< Dense > output) const
 
std::unique_ptr< Diagonal< ValueType > > extract_diagonal () const override
 Extracts the diagonal entries of the matrix into a vector. More...
 
void extract_diagonal (ptr_param< Diagonal< ValueType >> output) const
 Writes the diagonal of this matrix into an existing diagonal matrix. More...
 
std::unique_ptr< absolute_type > compute_absolute () const override
 Gets the AbsoluteLinOp. More...
 
void compute_absolute (ptr_param< absolute_type > output) const
 Writes the absolute values of this matrix into an existing matrix. More...
 
void compute_absolute_inplace () override
 Compute absolute inplace on each element.
 
std::unique_ptr< complex_type > make_complex () const
 Creates a complex copy of the original matrix. More...
 
void make_complex (ptr_param< complex_type > result) const
 Writes a complex copy of the original matrix to a given complex matrix. More...
 
std::unique_ptr< real_type > get_real () const
 Creates a new real matrix and extracts the real part of the original matrix into that.
 
void get_real (ptr_param< real_type > result) const
 Extracts the real part of the original matrix into a given real matrix.
 
std::unique_ptr< real_type > get_imag () const
 Creates a new real matrix and extracts the imaginary part of the original matrix into that.
 
void get_imag (ptr_param< real_type > result) const
 Extracts the imaginary part of the original matrix into a given real matrix.
 
value_type * get_values () noexcept
 Returns a pointer to the array of values of the matrix. More...
 
const value_type * get_const_values () const noexcept
 Returns a pointer to the array of values of the matrix. More...
 
size_type get_stride () const noexcept
 Returns the stride of the matrix. More...
 
size_type get_num_stored_elements () const noexcept
 Returns the number of elements explicitly stored in the matrix. More...
 
value_type & at (size_type row, size_type col) noexcept
 Returns a single element of the matrix. More...
 
value_type at (size_type row, size_type col) const noexcept
 Returns a single element of the matrix. More...
 
ValueType & at (size_type idx) noexcept
 Returns a single element of the matrix. More...
 
ValueType at (size_type idx) const noexcept
 Returns a single element of the matrix. More...
 
void scale (ptr_param< const LinOp > alpha)
 Scales the matrix with a scalar (aka: BLAS scal). More...
 
void inv_scale (ptr_param< const LinOp > alpha)
 Scales the matrix with the inverse of a scalar. More...
 
void add_scaled (ptr_param< const LinOp > alpha, ptr_param< const LinOp > b)
 Adds b scaled by alpha to the matrix (aka: BLAS axpy). More...
 
void sub_scaled (ptr_param< const LinOp > alpha, ptr_param< const LinOp > b)
 Subtracts b scaled by alpha from the matrix (aka: BLAS axpy). More...
 
void compute_dot (ptr_param< const LinOp > b, ptr_param< LinOp > result) const
 Computes the column-wise dot product of this matrix and b. More...
 
void compute_dot (ptr_param< const LinOp > b, ptr_param< LinOp > result, array< char > &tmp) const
 Computes the column-wise dot product of this matrix and b. More...
 
void compute_conj_dot (ptr_param< const LinOp > b, ptr_param< LinOp > result) const
 Computes the column-wise dot product of conj(this matrix) and b. More...
 
void compute_conj_dot (ptr_param< const LinOp > b, ptr_param< LinOp > result, array< char > &tmp) const
 Computes the column-wise dot product of conj(this matrix) and b. More...
 
void compute_norm2 (ptr_param< LinOp > result) const
 Computes the column-wise Euclidean (L^2) norm of this matrix. More...
 
void compute_norm2 (ptr_param< LinOp > result, array< char > &tmp) const
 Computes the column-wise Euclidean (L^2) norm of this matrix. More...
 
void compute_norm1 (ptr_param< LinOp > result) const
 Computes the column-wise (L^1) norm of this matrix. More...
 
void compute_norm1 (ptr_param< LinOp > result, array< char > &tmp) const
 Computes the column-wise (L^1) norm of this matrix. More...
 
void compute_squared_norm2 (ptr_param< LinOp > result) const
 Computes the square of the column-wise Euclidean (L^2) norm of this matrix. More...
 
void compute_squared_norm2 (ptr_param< LinOp > result, array< char > &tmp) const
 Computes the square of the column-wise Euclidean (L^2) norm of this matrix. More...
 
void compute_mean (ptr_param< LinOp > result) const
 Computes the column-wise arithmetic mean of this matrix. More...
 
void compute_mean (ptr_param< LinOp > result, array< char > &tmp) const
 Computes the column-wise arithmetic mean of this matrix. More...
 
std::unique_ptr< Densecreate_submatrix (const span &rows, const span &columns, const size_type stride)
 Create a submatrix from the original matrix. More...
 
std::unique_ptr< Densecreate_submatrix (const span &rows, const span &columns)
 Create a submatrix from the original matrix. More...
 
std::unique_ptr< real_type > create_real_view ()
 Create a real view of the (potentially) complex original matrix. More...
 
std::unique_ptr< const real_type > create_real_view () const
 Create a real view of the (potentially) complex original matrix. More...
 
Denseoperator= (const Dense &)
 Copy-assigns a Dense matrix. More...
 
Denseoperator= (Dense &&)
 Move-assigns a Dense matrix. More...
 
 Dense (const Dense &)
 Copy-constructs a Dense matrix. More...
 
 Dense (Dense &&)
 Move-constructs a Dense matrix. More...
 
- Public Member Functions inherited from gko::EnableLinOp< Dense< ValueType > >
const Dense< ValueType > * apply (ptr_param< const LinOp > b, ptr_param< LinOp > x) const
 
Dense< ValueType > * apply (ptr_param< const LinOp > b, ptr_param< LinOp > x)
 
const Dense< ValueType > * apply (ptr_param< const LinOp > alpha, ptr_param< const LinOp > b, ptr_param< const LinOp > beta, ptr_param< LinOp > x) const
 
Dense< ValueType > * 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< Dense< ValueType >, LinOp >
std::unique_ptr< Dense< ValueType > > create_default (std::shared_ptr< const Executor > exec) const
 
std::unique_ptr< Dense< ValueType > > create_default () const
 
std::unique_ptr< Dense< ValueType > > clone (std::shared_ptr< const Executor > exec) const
 
std::unique_ptr< Dense< ValueType > > clone () const
 
Dense< ValueType > * copy_from (const PolymorphicObject *other)
 
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, Dense< ValueType > > * copy_from (std::unique_ptr< Derived > &&other)
 
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, Dense< ValueType > > * copy_from (const std::unique_ptr< Derived > &other)
 
Dense< ValueType > * copy_from (const std::shared_ptr< const PolymorphicObject > &other)
 
Dense< ValueType > * move_from (ptr_param< PolymorphicObject > other)
 
Dense< ValueType > * 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< Dense< ValueType > >
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< Dense< ValueType > >
void convert_to (ptr_param< result_type > result) const
 
void move_to (ptr_param< result_type > result)
 
- Public Member Functions inherited from gko::ConvertibleTo< Dense< next_precision< 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, int32 > >
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, int64 > >
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< Csr< ValueType, int32 > >
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< Csr< ValueType, int64 > >
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, int32 > >
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, int64 > >
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, int32 > >
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, int64 > >
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, int32 > >
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, int64 > >
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, int32 > >
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, int64 > >
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, int32 > >
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, int64 > >
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, int32 >
virtual void read (const matrix_data< ValueType, int32 > &data)=0
 Reads a matrix from a matrix_data structure. More...
 
void read (const matrix_assembly_data< ValueType, int32 > &data)
 Reads a matrix from a matrix_assembly_data structure. More...
 
virtual void read (const device_matrix_data< ValueType, int32 > &data)
 Reads a matrix from a device_matrix_data structure. More...
 
virtual void read (device_matrix_data< ValueType, int32 > &&data)
 Reads a matrix from a device_matrix_data structure. More...
 
- Public Member Functions inherited from gko::ReadableFromMatrixData< ValueType, int64 >
virtual void read (const matrix_data< ValueType, int64 > &data)=0
 Reads a matrix from a matrix_data structure. More...
 
void read (const matrix_assembly_data< ValueType, int64 > &data)
 Reads a matrix from a matrix_assembly_data structure. More...
 
virtual void read (const device_matrix_data< ValueType, int64 > &data)
 Reads a matrix from a device_matrix_data structure. More...
 
virtual void read (device_matrix_data< ValueType, int64 > &&data)
 Reads a matrix from a device_matrix_data structure. More...
 
- Public Member Functions inherited from gko::WritableToMatrixData< ValueType, int32 >
virtual void write (matrix_data< ValueType, int32 > &data) const=0
 Writes a matrix to a matrix_data structure. More...
 
- Public Member Functions inherited from gko::WritableToMatrixData< ValueType, int64 >
virtual void write (matrix_data< ValueType, int64 > &data) const=0
 Writes a matrix to a matrix_data structure. More...
 
- Public Member Functions inherited from gko::EnableAbsoluteComputation< remove_complex< Dense< ValueType > > >
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< Densecreate_with_config_of (ptr_param< const Dense > other)
 Creates a Dense matrix with the same size and stride as another Dense matrix. More...
 
static std::unique_ptr< Densecreate_with_type_of (ptr_param< const Dense > other, std::shared_ptr< const Executor > exec, const dim< 2 > &size=dim< 2 >{})
 Creates a Dense matrix with the same type as another Dense matrix but on a different executor and with a different size. More...
 
static std::unique_ptr< Densecreate_with_type_of (ptr_param< const Dense > other, std::shared_ptr< const Executor > exec, const dim< 2 > &size, size_type stride)
 
static std::unique_ptr< Densecreate_with_type_of (ptr_param< const Dense > other, std::shared_ptr< const Executor > exec, const dim< 2 > &size, const dim< 2 > &local_size, size_type stride)
 
static std::unique_ptr< Densecreate_view_of (ptr_param< Dense > other)
 Creates a Dense matrix, where the underlying array is a view of another Dense matrix' array. More...
 
static std::unique_ptr< const Densecreate_const_view_of (ptr_param< const Dense > other)
 Creates a immutable Dense matrix, where the underlying array is a view of another Dense matrix' array. More...
 
static std::unique_ptr< const Densecreate_const (std::shared_ptr< const Executor > exec, const dim< 2 > &size, gko::detail::const_array_view< ValueType > &&values, size_type stride)
 Creates a constant (immutable) Dense matrix from a constant array. More...
 
- Static Public Member Functions inherited from gko::EnableCreateMethod< Dense< ValueType > >
static std::unique_ptr< Dense< ValueType > > create (Args &&... args)
 

Friends

class EnableCreateMethod< Dense >
 
class EnablePolymorphicObject< Dense, LinOp >
 
class Coo< ValueType, int32 >
 
class Coo< ValueType, int64 >
 
class Csr< ValueType, int32 >
 
class Csr< ValueType, int64 >
 
class Diagonal< ValueType >
 
class Ell< ValueType, int32 >
 
class Ell< ValueType, int64 >
 
class Fbcsr< ValueType, int32 >
 
class Fbcsr< ValueType, int64 >
 
class Hybrid< ValueType, int32 >
 
class Hybrid< ValueType, int64 >
 
class Sellp< ValueType, int32 >
 
class Sellp< ValueType, int64 >
 
class SparsityCsr< ValueType, int32 >
 
class SparsityCsr< ValueType, int64 >
 
class Dense< to_complex< ValueType > >
 
class experimental::distributed::Vector< ValueType >
 
class Dense< next_precision< ValueType > >
 

Detailed Description

template<typename ValueType = default_precision>
class gko::matrix::Dense< ValueType >

Dense is a matrix format which explicitly stores all values of the matrix.

The values are stored in row-major format (values belonging to the same row appear consecutive in the memory). Optionally, rows can be padded for better memory access.

Template Parameters
ValueTypeprecision of matrix elements
Note
While this format is not very useful for storing sparse matrices, it is often suitable to store vectors, and sets of vectors.

Constructor & Destructor Documentation

◆ Dense() [1/2]

template<typename ValueType = default_precision>
gko::matrix::Dense< ValueType >::Dense ( const Dense< ValueType > &  )

Copy-constructs a Dense matrix.

Inherits executor and dimensions, but copies data without padding.

◆ Dense() [2/2]

template<typename ValueType = default_precision>
gko::matrix::Dense< ValueType >::Dense ( Dense< ValueType > &&  )

Move-constructs a Dense matrix.

Inherits executor, dimensions and data with padding. The moved-from object is empty (0x0 with empty Array).

Member Function Documentation

◆ add_scaled()

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::add_scaled ( ptr_param< const LinOp alpha,
ptr_param< const LinOp b 
)

Adds b scaled by alpha to the matrix (aka: BLAS axpy).

Parameters
alphaIf alpha is 1x1 Dense matrix, the entire matrix is scaled by alpha. If it is a Dense row vector of values, then i-th column of the matrix is scaled with the i-th element of alpha (the number of columns of alpha has to match the number of columns of the matrix).
ba matrix of the same dimension as this

◆ at() [1/4]

template<typename ValueType = default_precision>
ValueType gko::matrix::Dense< ValueType >::at ( size_type  idx) const
inlinenoexcept

Returns a single element of the matrix.

Useful for iterating across all elements of the matrix. However, it is less efficient than the two-parameter variant of this method.

Parameters
idxa linear index of the requested element (ignoring the stride)
Note
the method has to be called on the same Executor the matrix is stored at (e.g. trying to call this method on a GPU matrix from the OMP results in a runtime error)

◆ at() [2/4]

template<typename ValueType = default_precision>
ValueType& gko::matrix::Dense< ValueType >::at ( size_type  idx)
inlinenoexcept

Returns a single element of the matrix.

Useful for iterating across all elements of the matrix. However, it is less efficient than the two-parameter variant of this method.

Parameters
idxa linear index of the requested element (ignoring the stride)
Note
the method has to be called on the same Executor the matrix is stored at (e.g. trying to call this method on a GPU matrix from the OMP results in a runtime error)

◆ at() [3/4]

template<typename ValueType = default_precision>
value_type gko::matrix::Dense< ValueType >::at ( size_type  row,
size_type  col 
) const
inlinenoexcept

Returns a single element of the matrix.

Parameters
rowthe row of the requested element
colthe column of the requested element
Note
the method has to be called on the same Executor the matrix is stored at (e.g. trying to call this method on a GPU matrix from the OMP results in a runtime error)

◆ at() [4/4]

template<typename ValueType = default_precision>
value_type& gko::matrix::Dense< ValueType >::at ( size_type  row,
size_type  col 
)
inlinenoexcept

Returns a single element of the matrix.

Parameters
rowthe row of the requested element
colthe column of the requested element
Note
the method has to be called on the same Executor the matrix is stored at (e.g. trying to call this method on a GPU matrix from the OMP results in a runtime error)

Referenced by gko::initialize().

◆ column_permute() [1/4]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::column_permute ( const array< int32 > *  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< int32 >.

◆ column_permute() [2/4]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::column_permute ( const array< int32 > *  permutation_indices,
ptr_param< Dense< ValueType > >  output 
) const

Writes the column-permuted matrix into the given output matrix.

Parameters
permutation_indicesThe array containing permutation indices. It must have this->get_size()[1] elements.
outputThe output matrix. It must have the dimensions this->get_size()
See also
Dense::column_permute(const array<int32>*)

◆ column_permute() [3/4]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::column_permute ( const array< int64 > *  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< int64 >.

◆ column_permute() [4/4]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::column_permute ( const array< int64 > *  permutation_indices,
ptr_param< Dense< ValueType > >  output 
) const

◆ compute_absolute() [1/2]

template<typename ValueType = default_precision>
std::unique_ptr<absolute_type> gko::matrix::Dense< ValueType >::compute_absolute ( ) const
overridevirtual

Gets the AbsoluteLinOp.

Returns
a pointer to the new absolute object

Implements gko::EnableAbsoluteComputation< remove_complex< Dense< ValueType > > >.

◆ compute_absolute() [2/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_absolute ( ptr_param< absolute_type >  output) const

Writes the absolute values of this matrix into an existing matrix.

Parameters
outputThe output matrix. Its size must match the size of this matrix.
See also
Dense::compute_absolute()

◆ compute_conj_dot() [1/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_conj_dot ( ptr_param< const LinOp b,
ptr_param< LinOp result 
) const

Computes the column-wise dot product of conj(this matrix) and b.

Parameters
ba Dense matrix of same dimension as this
resulta Dense row vector, used to store the dot product (the number of column in the vector must match the number of columns of this)

◆ compute_conj_dot() [2/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_conj_dot ( ptr_param< const LinOp b,
ptr_param< LinOp result,
array< char > &  tmp 
) const

Computes the column-wise dot product of conj(this matrix) and b.

Parameters
ba Dense matrix of same dimension as this
resulta Dense row vector, used to store the dot product (the number of column in the vector must match the number of columns of this)
tmpthe temporary storage to use for partial sums during the reduction computation. It may be resized and/or reset to the correct executor.

◆ compute_dot() [1/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_dot ( ptr_param< const LinOp b,
ptr_param< LinOp result 
) const

Computes the column-wise dot product of this matrix and b.

Parameters
ba Dense matrix of same dimension as this
resulta Dense row vector, used to store the dot product (the number of column in the vector must match the number of columns of this)

◆ compute_dot() [2/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_dot ( ptr_param< const LinOp b,
ptr_param< LinOp result,
array< char > &  tmp 
) const

Computes the column-wise dot product of this matrix and b.

Parameters
ba Dense matrix of same dimension as this
resulta Dense row vector, used to store the dot product (the number of column in the vector must match the number of columns of this)
tmpthe temporary storage to use for partial sums during the reduction computation. It may be resized and/or reset to the correct executor.

◆ compute_mean() [1/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_mean ( ptr_param< LinOp result) const

Computes the column-wise arithmetic mean of this matrix.

Parameters
resulta Dense row vector, used to store the mean (the number of columns in the vector must match the number of columns of this)

◆ compute_mean() [2/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_mean ( ptr_param< LinOp result,
array< char > &  tmp 
) const

Computes the column-wise arithmetic mean of this matrix.

Parameters
resulta Dense row vector, used to store the mean (the number of columns in the vector must match the number of columns of this)
tmpthe temporary storage to use for partial sums during the reduction computation. It may be resized and/or reset to the correct executor.

◆ compute_norm1() [1/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_norm1 ( ptr_param< LinOp result) const

Computes the column-wise (L^1) norm of this matrix.

Parameters
resulta Dense row vector, used to store the norm (the number of columns in the vector must match the number of columns of this)

◆ compute_norm1() [2/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_norm1 ( ptr_param< LinOp result,
array< char > &  tmp 
) const

Computes the column-wise (L^1) norm of this matrix.

Parameters
resulta Dense row vector, used to store the norm (the number of columns in the vector must match the number of columns of this)
tmpthe temporary storage to use for partial sums during the reduction computation. It may be resized and/or reset to the correct executor.

◆ compute_norm2() [1/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_norm2 ( ptr_param< LinOp result) const

Computes the column-wise Euclidean (L^2) norm of this matrix.

Parameters
resulta Dense row vector, used to store the norm (the number of columns in the vector must match the number of columns of this)

◆ compute_norm2() [2/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_norm2 ( ptr_param< LinOp result,
array< char > &  tmp 
) const

Computes the column-wise Euclidean (L^2) norm of this matrix.

Parameters
resulta Dense row vector, used to store the norm (the number of columns in the vector must match the number of columns of this)
tmpthe temporary storage to use for partial sums during the reduction computation. It may be resized and/or reset to the correct executor.

◆ compute_squared_norm2() [1/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_squared_norm2 ( ptr_param< LinOp result) const

Computes the square of the column-wise Euclidean (L^2) norm of this matrix.

Parameters
resulta Dense row vector, used to store the norm (the number of columns in the vector must match the number of columns of this)

◆ compute_squared_norm2() [2/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::compute_squared_norm2 ( ptr_param< LinOp result,
array< char > &  tmp 
) const

Computes the square of the column-wise Euclidean (L^2) norm of this matrix.

Parameters
resulta Dense row vector, used to store the norm (the number of columns in the vector must match the number of columns of this)
tmpthe temporary storage to use for partial sums during the reduction computation. It may be resized and/or reset to the correct executor.

◆ conj_transpose() [1/2]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::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.

◆ conj_transpose() [2/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::conj_transpose ( ptr_param< Dense< ValueType > >  output) const

Writes the conjugate-transposed matrix into the given output matrix.

Parameters
outputThe output matrix. It must have the dimensions gko::transpose(this->get_size())

◆ create_const()

template<typename ValueType = default_precision>
static std::unique_ptr<const Dense> gko::matrix::Dense< ValueType >::create_const ( std::shared_ptr< const Executor exec,
const dim< 2 > &  size,
gko::detail::const_array_view< ValueType > &&  values,
size_type  stride 
)
inlinestatic

Creates a constant (immutable) Dense matrix from a constant array.

Parameters
execthe executor to create the matrix on
sizethe dimensions of the matrix
valuesthe value array of the matrix
stridethe row-stride of the matrix
Returns
A smart pointer to the constant matrix wrapping the input array (if it resides on the same executor as the matrix) or a copy of the array on the correct executor.

◆ create_const_view_of()

template<typename ValueType = default_precision>
static std::unique_ptr<const Dense> gko::matrix::Dense< ValueType >::create_const_view_of ( ptr_param< const Dense< ValueType > >  other)
inlinestatic

Creates a immutable Dense matrix, where the underlying array is a view of another Dense matrix' array.

Parameters
otherThe other matrix on which to create the view
Returns
A immutable Dense matrix that is a view of other

Referenced by gko::make_const_dense_view().

◆ create_real_view() [1/2]

template<typename ValueType = default_precision>
std::unique_ptr<real_type> gko::matrix::Dense< ValueType >::create_real_view ( )

Create a real view of the (potentially) complex original matrix.

If the original matrix is real, nothing changes. If the original matrix is complex, the result is created by viewing the complex matrix with as real with a reinterpret_cast with twice the number of columns and double the stride.

◆ create_real_view() [2/2]

template<typename ValueType = default_precision>
std::unique_ptr<const real_type> gko::matrix::Dense< ValueType >::create_real_view ( ) const

Create a real view of the (potentially) complex original matrix.

If the original matrix is real, nothing changes. If the original matrix is complex, the result is created by viewing the complex matrix with as real with a reinterpret_cast with twice the number of columns and double the stride.

◆ create_submatrix() [1/2]

template<typename ValueType = default_precision>
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::create_submatrix ( const span rows,
const span columns 
)
inline

Create a submatrix from the original matrix.

Parameters
rowsrow span
columnscolumn span

◆ create_submatrix() [2/2]

template<typename ValueType = default_precision>
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::create_submatrix ( const span rows,
const span columns,
const size_type  stride 
)
inline

Create a submatrix from the original matrix.

Warning: defining stride for this create_submatrix method might cause wrong memory access. Better use the create_submatrix(rows, columns) method instead.

Parameters
rowsrow span
columnscolumn span
stridestride of the new submatrix.

Referenced by gko::matrix::Dense< value_type >::create_submatrix().

◆ create_view_of()

template<typename ValueType = default_precision>
static std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::create_view_of ( ptr_param< Dense< ValueType > >  other)
inlinestatic

Creates a Dense matrix, where the underlying array is a view of another Dense matrix' array.

Parameters
otherThe other matrix on which to create the view
Returns
A Dense matrix that is a view of other

Referenced by gko::make_dense_view().

◆ create_with_config_of()

template<typename ValueType = default_precision>
static std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::create_with_config_of ( ptr_param< const Dense< ValueType > >  other)
inlinestatic

Creates a Dense matrix with the same size and stride as another Dense matrix.

Parameters
otherThe other matrix whose configuration needs to copied.

◆ create_with_type_of() [1/3]

template<typename ValueType = default_precision>
static std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::create_with_type_of ( ptr_param< const Dense< ValueType > >  other,
std::shared_ptr< const Executor exec,
const dim< 2 > &  size,
const dim< 2 > &  local_size,
size_type  stride 
)
inlinestatic

Parameters
local_sizeUnused
strideThe stride of the new matrix.
Note
This is an overload to stay consistent with gko::experimental::distributed::Vector

◆ create_with_type_of() [2/3]

template<typename ValueType = default_precision>
static std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::create_with_type_of ( ptr_param< const Dense< ValueType > >  other,
std::shared_ptr< const Executor exec,
const dim< 2 > &  size,
size_type  stride 
)
inlinestatic

Parameters
strideThe stride of the new matrix.
Note
This is an overload which allows full parameter specification.

◆ create_with_type_of() [3/3]

template<typename ValueType = default_precision>
static std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::create_with_type_of ( ptr_param< const Dense< ValueType > >  other,
std::shared_ptr< const Executor exec,
const dim< 2 > &  size = dim<2>{} 
)
inlinestatic

Creates a Dense matrix with the same type as another Dense matrix but on a different executor and with a different size.

Parameters
otherThe other matrix whose type we target.
execThe executor of the new matrix.
sizeThe size of the new matrix.
strideThe stride of the new matrix.
Returns
a Dense matrix with the type of other.

◆ extract_diagonal() [1/2]

template<typename ValueType = default_precision>
std::unique_ptr<Diagonal<ValueType> > gko::matrix::Dense< ValueType >::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 >.

◆ extract_diagonal() [2/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::extract_diagonal ( ptr_param< Diagonal< ValueType >>  output) const

Writes the diagonal of this matrix into an existing diagonal matrix.

Parameters
outputThe output matrix. Its size must match the size of this matrix's diagonal.
See also
Dense::extract_diagonal()

◆ fill()

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::fill ( const ValueType  value)

Fill the dense matrix with a given value.

Parameters
valuethe value to be filled

◆ get_const_values()

template<typename ValueType = default_precision>
const value_type* gko::matrix::Dense< ValueType >::get_const_values ( ) const
inlinenoexcept

Returns a pointer to the array of values of the matrix.

Returns
the pointer to the array of values
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.

◆ get_num_stored_elements()

template<typename ValueType = default_precision>
size_type gko::matrix::Dense< ValueType >::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

◆ get_stride()

template<typename ValueType = default_precision>
size_type gko::matrix::Dense< ValueType >::get_stride ( ) const
inlinenoexcept

Returns the stride of the matrix.

Returns
the stride of the matrix.

Referenced by gko::matrix::Dense< value_type >::create_submatrix().

◆ get_values()

template<typename ValueType = default_precision>
value_type* gko::matrix::Dense< ValueType >::get_values ( )
inlinenoexcept

Returns a pointer to the array of values of the matrix.

Returns
the pointer to the array of values

◆ inv_scale()

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::inv_scale ( ptr_param< const LinOp alpha)

Scales the matrix with the inverse of a scalar.

Parameters
alphaIf alpha is 1x1 Dense matrix, the entire matrix is scaled by 1 / alpha. If it is a Dense row vector of values, then i-th column of the matrix is scaled with the inverse of the i-th element of alpha (the number of columns of alpha has to match the number of columns of the matrix).

◆ inverse_column_permute() [1/4]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::inverse_column_permute ( const array< int32 > *  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< int32 >.

◆ inverse_column_permute() [2/4]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::inverse_column_permute ( const array< int32 > *  permutation_indices,
ptr_param< Dense< ValueType > >  output 
) const

Writes the inverse column-permuted matrix into the given output matrix.

Parameters
permutation_indicesThe array containing permutation indices. It must have this->get_size()[1] elements.
outputThe output matrix. It must have the dimensions this->get_size()
See also
Dense::inverse_column_permute(const array<int32>*)

◆ inverse_column_permute() [3/4]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::inverse_column_permute ( const array< int64 > *  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< int64 >.

◆ inverse_column_permute() [4/4]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::inverse_column_permute ( const array< int64 > *  permutation_indices,
ptr_param< Dense< ValueType > >  output 
) const

◆ inverse_permute() [1/4]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::inverse_permute ( const array< int32 > *  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< int32 >.

◆ inverse_permute() [2/4]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::inverse_permute ( const array< int32 > *  permutation_indices,
ptr_param< Dense< ValueType > >  output 
) const

Writes the inverse symmetrically permuted matrix into the given output matrix.

Parameters
permutation_indicesThe array containing permutation indices. It must have this->get_size()[0] elements.
outputThe output matrix. It must have the dimensions this->get_size()
See also
Dense::inverse_permute(const array<int32>*)

◆ inverse_permute() [3/4]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::inverse_permute ( const array< int64 > *  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< int64 >.

◆ inverse_permute() [4/4]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::inverse_permute ( const array< int64 > *  permutation_indices,
ptr_param< Dense< ValueType > >  output 
) const

◆ inverse_row_permute() [1/4]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::inverse_row_permute ( const array< int32 > *  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< int32 >.

◆ inverse_row_permute() [2/4]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::inverse_row_permute ( const array< int32 > *  permutation_indices,
ptr_param< Dense< ValueType > >  output 
) const

Writes the inverse row-permuted matrix into the given output matrix.

Parameters
permutation_indicesThe array containing permutation indices. It must have this->get_size()[0] elements.
outputThe output matrix. It must have the dimensions this->get_size()
See also
Dense::inverse_row_permute(const array<int32>*)

◆ inverse_row_permute() [3/4]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::inverse_row_permute ( const array< int64 > *  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< int64 >.

◆ inverse_row_permute() [4/4]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::inverse_row_permute ( const array< int64 > *  permutation_indices,
ptr_param< Dense< ValueType > >  output 
) const

◆ make_complex() [1/2]

template<typename ValueType = default_precision>
std::unique_ptr<complex_type> gko::matrix::Dense< ValueType >::make_complex ( ) const

Creates a complex copy of the original matrix.

If the original matrix was real, the imaginary part of the result will be zero.

◆ make_complex() [2/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::make_complex ( ptr_param< complex_type >  result) const

Writes a complex copy of the original matrix to a given complex matrix.

If the original matrix was real, the imaginary part of the result will be zero.

◆ operator=() [1/2]

template<typename ValueType = default_precision>
Dense& gko::matrix::Dense< ValueType >::operator= ( const Dense< ValueType > &  )

Copy-assigns a Dense matrix.

Preserves the executor, reallocates the matrix with minimal stride if the dimensions don't match, then copies the data over, ignoring padding.

◆ operator=() [2/2]

template<typename ValueType = default_precision>
Dense& gko::matrix::Dense< ValueType >::operator= ( Dense< ValueType > &&  )

Move-assigns a Dense 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).

◆ permute() [1/12]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::permute ( const array< int32 > *  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< int32 >.

◆ permute() [2/12]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::permute ( const array< int32 > *  permutation_indices,
ptr_param< Dense< ValueType > >  output 
) const

Writes the symmetrically permuted matrix into the given output matrix.

Parameters
permutation_indicesThe array containing permutation indices. It must have this->get_size()[0] elements.
outputThe output matrix. It must have the dimensions this->get_size()
See also
Dense::permute(const array<int32>*)

◆ permute() [3/12]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::permute ( const array< int64 > *  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< int64 >.

◆ permute() [4/12]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::permute ( const array< int64 > *  permutation_indices,
ptr_param< Dense< ValueType > >  output 
) const

◆ permute() [5/12]

template<typename ValueType = default_precision>
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::permute ( ptr_param< const Permutation< int32 >>  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() [6/12]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::permute ( ptr_param< const Permutation< int32 >>  permutation,
ptr_param< Dense< ValueType > >  output,
permute_mode  mode 
) const

Overload of permute(ptr_param<const Permutation<int32>>, permute_mode) that writes the permuted copy into an existing Dense matrix.

Parameters
outputthe output matrix.

◆ permute() [7/12]

template<typename ValueType = default_precision>
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::permute ( ptr_param< const Permutation< int32 >>  row_permutation,
ptr_param< const Permutation< int32 >>  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() [8/12]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::permute ( ptr_param< const Permutation< int32 >>  row_permutation,
ptr_param< const Permutation< int32 >>  column_permutation,
ptr_param< Dense< ValueType > >  output,
bool  invert = false 
) const

Overload of permute(ptr_param<const Permutation<int32>>, ptr_param<const Permutation<int32>>, permute_mode) that writes the permuted copy into an existing Dense matrix.

Parameters
outputthe output matrix.

◆ permute() [9/12]

template<typename ValueType = default_precision>
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::permute ( ptr_param< const Permutation< int64 >>  permutation,
permute_mode  mode = permute_mode::symmetric 
) const

◆ permute() [10/12]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::permute ( ptr_param< const Permutation< int64 >>  permutation,
ptr_param< Dense< ValueType > >  output,
permute_mode  mode 
) const

◆ permute() [11/12]

template<typename ValueType = default_precision>
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::permute ( ptr_param< const Permutation< int64 >>  row_permutation,
ptr_param< const Permutation< int64 >>  column_permutation,
bool  invert = false 
) const

◆ permute() [12/12]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::permute ( ptr_param< const Permutation< int64 >>  row_permutation,
ptr_param< const Permutation< int64 >>  column_permutation,
ptr_param< Dense< ValueType > >  output,
bool  invert = false 
) const

◆ row_gather() [1/6]

template<typename ValueType = default_precision>
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::row_gather ( const array< int32 > *  gather_indices) const

Create a Dense matrix consisting of the given rows from this matrix.

Parameters
gather_indicespointer to an array containing row indices from this matrix. It may contain duplicates.
Returns
Dense matrix on the same executor with the same number of columns and gather_indices->get_num_elems() rows containing the gathered rows from this matrix: output(i,j) = input(gather_indices(i), j)

◆ row_gather() [2/6]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::row_gather ( const array< int32 > *  gather_indices,
ptr_param< LinOp row_collection 
) const

Copies the given rows from this matrix into row_collection

Parameters
gather_indicespointer to an array containing row indices from this matrix. It may contain duplicates.
row_collectionpointer to a LinOp that will store the gathered rows: row_collection(i,j) = input(gather_indices(i), j) It must have the same number of columns as this matrix and gather_indices->get_num_elems() rows.

◆ row_gather() [3/6]

template<typename ValueType = default_precision>
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::row_gather ( const array< int64 > *  gather_indices) const

Create a Dense matrix consisting of the given rows from this matrix.

Parameters
gather_indicespointer to an array containing row indices from this matrix. It may contain duplicates.
Returns
Dense matrix on the same executor with the same number of columns and gather_indices->get_num_elems() rows containing the gathered rows from this matrix: output(i,j) = input(gather_indices(i), j)

◆ row_gather() [4/6]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::row_gather ( const array< int64 > *  gather_indices,
ptr_param< LinOp row_collection 
) const

◆ row_gather() [5/6]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::row_gather ( ptr_param< const LinOp alpha,
const array< int32 > *  gather_indices,
ptr_param< const LinOp beta,
ptr_param< LinOp row_collection 
) const

Copies the given rows from this matrix into row_collection with scaling.

Parameters
alphascaling the result of row gathering
gather_indicespointer to an array containing row indices from this matrix. It may contain duplicates.
betascaling the input row_collection
row_collectionpointer to a LinOp that will store the gathered rows: row_collection(i,j) = input(gather_indices(i), j) It must have the same number of columns as this matrix and gather_indices->get_num_elems() rows.

◆ row_gather() [6/6]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::row_gather ( ptr_param< const LinOp alpha,
const array< int64 > *  gather_indices,
ptr_param< const LinOp beta,
ptr_param< LinOp row_collection 
) const

◆ row_permute() [1/4]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::row_permute ( const array< int32 > *  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< int32 >.

◆ row_permute() [2/4]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::row_permute ( const array< int32 > *  permutation_indices,
ptr_param< Dense< ValueType > >  output 
) const

Writes the row-permuted matrix into the given output matrix.

Parameters
permutation_indicesThe array containing permutation indices. It must have this->get_size()[0] elements.
outputThe output matrix. It must have the dimensions this->get_size()
See also
Dense::row_permute(const array<int32>*)

◆ row_permute() [3/4]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::row_permute ( const array< int64 > *  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< int64 >.

◆ row_permute() [4/4]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::row_permute ( const array< int64 > *  permutation_indices,
ptr_param< Dense< ValueType > >  output 
) const

◆ scale()

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::scale ( ptr_param< const LinOp alpha)

Scales the matrix with a scalar (aka: BLAS scal).

Parameters
alphaIf alpha is 1x1 Dense matrix, the entire matrix is scaled by alpha. If it is a Dense row vector of values, then i-th column of the matrix is scaled with the i-th element of alpha (the number of columns of alpha has to match the number of columns of the matrix).

◆ scale_permute() [1/8]

template<typename ValueType = default_precision>
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::scale_permute ( ptr_param< const ScaledPermutation< value_type, int32 >>  permutation,
permute_mode  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/8]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::scale_permute ( ptr_param< const ScaledPermutation< value_type, int32 >>  permutation,
ptr_param< Dense< ValueType > >  output,
permute_mode  mode 
) const

Overload of scale_permute(ptr_param<const ScaledPermutation<value_type, int32>>, permute_mode) that writes the permuted copy into an existing Dense matrix.

Parameters
outputthe output matrix.

◆ scale_permute() [3/8]

template<typename ValueType = default_precision>
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::scale_permute ( ptr_param< const ScaledPermutation< value_type, int32 >>  row_permutation,
ptr_param< const ScaledPermutation< value_type, int32 >>  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.

◆ scale_permute() [4/8]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::scale_permute ( ptr_param< const ScaledPermutation< value_type, int32 >>  row_permutation,
ptr_param< const ScaledPermutation< value_type, int32 >>  column_permutation,
ptr_param< Dense< ValueType > >  output,
bool  invert = false 
) const

Overload of scale_permute(ptr_param<const ScaledPermutation<value_type, int32>>, ptr_param<const ScaledPermutation<value_type, int32>>, bool) that writes the permuted copy into an existing Dense matrix.

Parameters
outputthe output matrix.

◆ scale_permute() [5/8]

template<typename ValueType = default_precision>
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::scale_permute ( ptr_param< const ScaledPermutation< value_type, int64 >>  permutation,
permute_mode  mode = permute_mode::symmetric 
) const

◆ scale_permute() [6/8]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::scale_permute ( ptr_param< const ScaledPermutation< value_type, int64 >>  permutation,
ptr_param< Dense< ValueType > >  output,
permute_mode  mode 
) const

◆ scale_permute() [7/8]

template<typename ValueType = default_precision>
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::scale_permute ( ptr_param< const ScaledPermutation< value_type, int64 >>  row_permutation,
ptr_param< const ScaledPermutation< value_type, int64 >>  column_permutation,
bool  invert = false 
) const

◆ scale_permute() [8/8]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::scale_permute ( ptr_param< const ScaledPermutation< value_type, int64 >>  row_permutation,
ptr_param< const ScaledPermutation< value_type, int64 >>  column_permutation,
ptr_param< Dense< ValueType > >  output,
bool  invert = false 
) const

◆ sub_scaled()

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::sub_scaled ( ptr_param< const LinOp alpha,
ptr_param< const LinOp b 
)

Subtracts b scaled by alpha from the matrix (aka: BLAS axpy).

Parameters
alphaIf alpha is 1x1 Dense matrix, b is scaled by alpha. If it is a Dense row vector of values, then i-th column of b is scaled with the i-th element of alpha (the number of columns of alpha has to match the number of columns of the matrix).
ba matrix of the same dimension as this

◆ transpose() [1/2]

template<typename ValueType = default_precision>
std::unique_ptr<LinOp> gko::matrix::Dense< ValueType >::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() [2/2]

template<typename ValueType = default_precision>
void gko::matrix::Dense< ValueType >::transpose ( ptr_param< Dense< ValueType > >  output) const

Writes the transposed matrix into the given output matrix.

Parameters
outputThe output matrix. It must have the dimensions gko::transpose(this->get_size())

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