![]() |
Ginkgo
Generated from pipelines/2118098289 branch based on develop. Ginkgo version 1.11.0
A numerical linear algebra library targeting many-core architectures
|
Dense is a matrix format which explicitly stores all values of the matrix. More...
#include <ginkgo/core/matrix/dense.hpp>
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 > > |
![]() | |
using | result_type = Dense< ValueType > |
![]() | |
using | result_type = Dense< ValueType > |
![]() | |
using | result_type = Dense< next_precision< ValueType > > |
![]() | |
using | result_type = Coo< ValueType, int32 > |
![]() | |
using | result_type = Coo< ValueType, int64 > |
![]() | |
using | result_type = Csr< ValueType, int32 > |
![]() | |
using | result_type = Csr< ValueType, int64 > |
![]() | |
using | result_type = Ell< ValueType, int32 > |
![]() | |
using | result_type = Ell< ValueType, int64 > |
![]() | |
using | result_type = Fbcsr< ValueType, int32 > |
![]() | |
using | result_type = Fbcsr< ValueType, int64 > |
![]() | |
using | result_type = Hybrid< ValueType, int32 > |
![]() | |
using | result_type = Hybrid< ValueType, int64 > |
![]() | |
using | result_type = Sellp< ValueType, int32 > |
![]() | |
using | result_type = Sellp< ValueType, int64 > |
![]() | |
using | result_type = SparsityCsr< ValueType, int32 > |
![]() | |
using | result_type = SparsityCsr< ValueType, int64 > |
![]() | |
using | value_type = ValueType |
![]() | |
using | value_type = ValueType |
using | index_type = int32 |
![]() | |
using | value_type = ValueType |
using | index_type = int64 |
![]() | |
using | value_type = ValueType |
using | index_type = int32 |
![]() | |
using | value_type = ValueType |
using | index_type = int64 |
![]() | |
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< LinOp > | transpose () const override |
Returns a LinOp representing the transpose of the Transposable object. More... | |
std::unique_ptr< LinOp > | conj_transpose () const override |
Returns a LinOp representing the conjugate transpose of the Transposable object. More... | |
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< Dense > | permute (ptr_param< const Permutation< int32 >> permutation, permute_mode mode=permute_mode::symmetric) const |
Creates a permuted copy ![]() ![]() ![]() | |
std::unique_ptr< Dense > | permute (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< Dense > | 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 ![]() ![]() ![]() ![]() | |
std::unique_ptr< Dense > | permute (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< Dense > | 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. More... | |
std::unique_ptr< Dense > | scale_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< Dense > | 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. More... | |
std::unique_ptr< Dense > | scale_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< LinOp > | permute (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< LinOp > | permute (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< LinOp > | inverse_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< LinOp > | inverse_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< LinOp > | row_permute (const array< int32 > *permutation_indices) const override |
Returns a LinOp representing the row permutation of the Permutable object. More... | |
std::unique_ptr< LinOp > | row_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< Dense > | row_gather (const array< int32 > *gather_indices) const |
Create a Dense matrix consisting of the given rows from this matrix. More... | |
std::unique_ptr< Dense > | row_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< LinOp > | column_permute (const array< int32 > *permutation_indices) const override |
Returns a LinOp representing the column permutation of the Permutable object. More... | |
std::unique_ptr< LinOp > | column_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< LinOp > | inverse_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< LinOp > | inverse_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< LinOp > | inverse_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< LinOp > | inverse_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< Dense > | create_submatrix (const span &rows, const span &columns, const size_type stride) |
Create a submatrix from the original matrix. More... | |
std::unique_ptr< Dense > | create_submatrix (const span &rows, const span &columns) |
Create a submatrix from the original matrix. More... | |
std::unique_ptr< Dense > | create_submatrix (const local_span &rows, const local_span &columns, dim< 2 > size) |
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... | |
Dense & | operator= (const Dense &) |
Copy-assigns a Dense matrix. More... | |
Dense & | operator= (Dense &&) |
Move-assigns a Dense matrix. More... | |
Dense (const Dense &) | |
Copy-constructs a Dense matrix. More... | |
Dense (Dense &&) | |
Move-constructs a Dense matrix. More... | |
![]() | |
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) |
![]() | |
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 () |
![]() | |
LinOp * | apply (ptr_param< const LinOp > b, ptr_param< LinOp > x) |
Applies a linear operator to a vector (or a sequence of vectors). More... | |
const LinOp * | apply (ptr_param< const LinOp > b, ptr_param< LinOp > x) const |
LinOp * | apply (ptr_param< const LinOp > alpha, ptr_param< const LinOp > b, ptr_param< const LinOp > beta, ptr_param< LinOp > x) |
Performs the operation x = alpha * op(b) + beta * x. More... | |
const LinOp * | apply (ptr_param< const LinOp > alpha, ptr_param< const LinOp > b, ptr_param< const LinOp > beta, ptr_param< LinOp > x) const |
const dim< 2 > & | get_size () const noexcept |
Returns the size of the operator. More... | |
virtual bool | apply_uses_initial_guess () const |
Returns true if the linear operator uses the data given in x as an initial guess. More... | |
LinOp & | operator= (const LinOp &)=default |
Copy-assigns a LinOp. More... | |
LinOp & | operator= (LinOp &&other) |
Move-assigns a LinOp. More... | |
LinOp (const LinOp &)=default | |
Copy-constructs a LinOp. More... | |
LinOp (LinOp &&other) | |
Move-constructs a LinOp. More... | |
![]() | |
std::unique_ptr< LinOp > | create_default (std::shared_ptr< const Executor > exec) const |
std::unique_ptr< LinOp > | create_default () const |
std::unique_ptr< LinOp > | clone (std::shared_ptr< const Executor > exec) const |
std::unique_ptr< LinOp > | clone () const |
LinOp * | copy_from (const PolymorphicObject *other) |
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, LinOp > * | copy_from (std::unique_ptr< Derived > &&other) |
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, LinOp > * | copy_from (const std::unique_ptr< Derived > &other) |
LinOp * | copy_from (const std::shared_ptr< const PolymorphicObject > &other) |
LinOp * | move_from (ptr_param< PolymorphicObject > other) |
LinOp * | clear () |
![]() | |
PolymorphicObject & | operator= (const PolymorphicObject &) |
std::unique_ptr< PolymorphicObject > | create_default (std::shared_ptr< const Executor > exec) const |
Creates a new "default" object of the same dynamic type as this object. More... | |
std::unique_ptr< PolymorphicObject > | create_default () const |
Creates a new "default" object of the same dynamic type as this object. More... | |
std::unique_ptr< PolymorphicObject > | clone (std::shared_ptr< const Executor > exec) const |
Creates a clone of the object. More... | |
std::unique_ptr< PolymorphicObject > | clone () const |
Creates a clone of the object. More... | |
PolymorphicObject * | copy_from (const PolymorphicObject *other) |
Copies another object into this object. More... | |
template<typename Derived , typename Deleter > | |
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, PolymorphicObject > * | copy_from (std::unique_ptr< Derived, Deleter > &&other) |
Moves another object into this object. More... | |
template<typename Derived , typename Deleter > | |
std::enable_if_t< std::is_base_of< PolymorphicObject, std::decay_t< Derived > >::value, PolymorphicObject > * | copy_from (const std::unique_ptr< Derived, Deleter > &other) |
Copies another object into this object. More... | |
PolymorphicObject * | copy_from (const std::shared_ptr< const PolymorphicObject > &other) |
Copies another object into this object. More... | |
PolymorphicObject * | move_from (ptr_param< PolymorphicObject > other) |
Moves another object into this object. More... | |
PolymorphicObject * | clear () |
Transforms the object into its default state. More... | |
std::shared_ptr< const Executor > | get_executor () const noexcept |
Returns the Executor of the object. More... | |
![]() | |
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 |
![]() | |
void | remove_logger (ptr_param< const Logger > logger) |
![]() | |
void | convert_to (result_type *result) const override |
Converts the implementer to an object of type result_type. More... | |
void | move_to (result_type *result) override |
Converts the implementer to an object of type result_type by moving data from this object. More... | |
![]() | |
void | convert_to (ptr_param< result_type > result) const |
void | move_to (ptr_param< result_type > result) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
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) |
![]() | |
std::unique_ptr< LinOp > | extract_diagonal_linop () const override |
Extracts the diagonal entries of the matrix into a vector. More... | |
![]() | |
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... | |
![]() | |
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... | |
![]() | |
virtual void | write (matrix_data< ValueType, int32 > &data) const=0 |
Writes a matrix to a matrix_data structure. More... | |
![]() | |
virtual void | write (matrix_data< ValueType, int64 > &data) const=0 |
Writes a matrix to a matrix_data structure. More... | |
![]() | |
std::unique_ptr< LinOp > | compute_absolute_linop () const override |
Gets the absolute LinOp. More... | |
![]() | |
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< Dense > | create_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< Dense > | create_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< Dense > | create_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< Dense > | create_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< Dense > | create_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 Dense > | create_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< Dense > | create (std::shared_ptr< const Executor > exec, const dim< 2 > &size={}, size_type stride=0) |
Creates an uninitialized Dense matrix of the specified size. More... | |
static std::unique_ptr< Dense > | create (std::shared_ptr< const Executor > exec, const dim< 2 > &size, array< value_type > values, size_type stride) |
Creates a Dense matrix from an already allocated (and initialized) array. More... | |
template<typename InputValueType > | |
static std::unique_ptr< Dense > | create (std::shared_ptr< const Executor > exec, const dim< 2 > &size, std::initializer_list< InputValueType > values, size_type stride) |
create(std::shared_ptr<const Executor>, More... | |
static std::unique_ptr< const Dense > | create_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... | |
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.
ValueType | precision of matrix elements |
gko::matrix::Dense< ValueType >::Dense | ( | const Dense< ValueType > & | ) |
Copy-constructs a Dense matrix.
Inherits executor and dimensions, but copies data without padding.
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).
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).
alpha | If 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). |
b | a matrix of the same dimension as this |
|
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.
idx | a linear index of the requested element (ignoring the stride) |
|
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.
idx | a linear index of the requested element (ignoring the stride) |
|
inlinenoexcept |
Returns a single element of the matrix.
row | the row of the requested element |
col | the column of the requested element |
|
inlinenoexcept |
Returns a single element of the matrix.
row | the row of the requested element |
col | the column of the requested element |
Referenced by gko::initialize().
|
overridevirtual |
Returns a LinOp representing the column permutation of the Permutable object.
In the resulting LinOp, the column i
contains the input column perm[i]
.
From the linear algebra perspective, with , this represents the operation
.
permutation_indices | the array of indices containing the permutation order perm . |
Implements gko::Permutable< int32 >.
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.
permutation_indices | The array containing permutation indices. It must have this->get_size()[1] elements. |
output | The output matrix. It must have the dimensions this->get_size() |
|
overridevirtual |
Returns a LinOp representing the column permutation of the Permutable object.
In the resulting LinOp, the column i
contains the input column perm[i]
.
From the linear algebra perspective, with , this represents the operation
.
permutation_indices | the array of indices containing the permutation order perm . |
Implements gko::Permutable< int64 >.
void gko::matrix::Dense< ValueType >::column_permute | ( | const array< int64 > * | permutation_indices, |
ptr_param< Dense< ValueType > > | output | ||
) | const |
|
overridevirtual |
Gets the AbsoluteLinOp.
Implements gko::EnableAbsoluteComputation< remove_complex< Dense< ValueType > > >.
void gko::matrix::Dense< ValueType >::compute_absolute | ( | ptr_param< absolute_type > | output | ) | const |
Writes the absolute values of this matrix into an existing matrix.
output | The output matrix. Its size must match the size of this matrix. |
void gko::matrix::Dense< ValueType >::compute_conj_dot | ( | ptr_param< const LinOp > | b, |
ptr_param< LinOp > | result | ||
) | const |
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
.
b | a Dense matrix of same dimension as this |
result | a Dense row vector, used to store the dot product (the number of column in the vector must match the number of columns of this) |
tmp | the temporary storage to use for partial sums during the reduction computation. It may be resized and/or reset to the correct executor. |
void gko::matrix::Dense< ValueType >::compute_dot | ( | ptr_param< const LinOp > | b, |
ptr_param< LinOp > | result | ||
) | const |
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
.
b | a Dense matrix of same dimension as this |
result | a Dense row vector, used to store the dot product (the number of column in the vector must match the number of columns of this) |
tmp | the temporary storage to use for partial sums during the reduction computation. It may be resized and/or reset to the correct executor. |
void gko::matrix::Dense< ValueType >::compute_mean | ( | ptr_param< LinOp > | result | ) | const |
Computes the column-wise arithmetic mean of this matrix.
result | a Dense row vector, used to store the mean (the number of columns in the vector must match the number of columns of this) |
void gko::matrix::Dense< ValueType >::compute_mean | ( | ptr_param< LinOp > | result, |
array< char > & | tmp | ||
) | const |
Computes the column-wise arithmetic mean of this matrix.
result | a Dense row vector, used to store the mean (the number of columns in the vector must match the number of columns of this) |
tmp | the temporary storage to use for partial sums during the reduction computation. It may be resized and/or reset to the correct executor. |
void gko::matrix::Dense< ValueType >::compute_norm1 | ( | ptr_param< LinOp > | result | ) | const |
Computes the column-wise (L^1) norm of this matrix.
result | a Dense row vector, used to store the norm (the number of columns in the vector must match the number of columns of this) |
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.
result | a Dense row vector, used to store the norm (the number of columns in the vector must match the number of columns of this) |
tmp | the temporary storage to use for partial sums during the reduction computation. It may be resized and/or reset to the correct executor. |
void gko::matrix::Dense< ValueType >::compute_norm2 | ( | ptr_param< LinOp > | result | ) | const |
Computes the column-wise Euclidean (L^2) norm of this matrix.
result | a Dense row vector, used to store the norm (the number of columns in the vector must match the number of columns of this) |
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.
result | a Dense row vector, used to store the norm (the number of columns in the vector must match the number of columns of this) |
tmp | the temporary storage to use for partial sums during the reduction computation. It may be resized and/or reset to the correct executor. |
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.
result | a Dense row vector, used to store the norm (the number of columns in the vector must match the number of columns of this) |
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.
result | a Dense row vector, used to store the norm (the number of columns in the vector must match the number of columns of this) |
tmp | the temporary storage to use for partial sums during the reduction computation. It may be resized and/or reset to the correct executor. |
|
overridevirtual |
Returns a LinOp representing the conjugate transpose of the Transposable object.
Implements gko::Transposable.
void gko::matrix::Dense< ValueType >::conj_transpose | ( | ptr_param< Dense< ValueType > > | output | ) | const |
Writes the conjugate-transposed matrix into the given output matrix.
output | The output matrix. It must have the dimensions gko::transpose(this->get_size()) |
|
static |
Creates a Dense matrix from an already allocated (and initialized) array.
exec | Executor associated to the matrix |
size | size of the matrix |
values | array of matrix values |
stride | stride of the rows (i.e. offset between the first elements of two consecutive rows, expressed as the number of matrix elements) |
values
is not an rvalue, not an array of ValueType, or is on the wrong executor, an internal copy will be created, and the original array data will not be used in the matrix.
|
inlinestatic |
create(std::shared_ptr<const Executor>,
create(std::shared_ptr<const Executor>, const dim<2>&, array<value_type>, size_type)
|
static |
Creates an uninitialized Dense matrix of the specified size.
exec | Executor associated to the matrix |
size | size of the matrix |
stride | stride of the rows (i.e. offset between the first elements of two consecutive rows, expressed as the number of matrix elements). If it is set to 0, size[1] will be used instead. |
Referenced by gko::matrix::Dense< value_type >::create().
|
static |
Creates a constant (immutable) Dense matrix from a constant array.
exec | the executor to create the matrix on |
size | the dimensions of the matrix |
values | the value array of the matrix |
stride | the row-stride of the matrix |
|
inlinestatic |
Creates a immutable Dense matrix, where the underlying array is a view of another Dense matrix' array.
other | The other matrix on which to create the view |
Referenced by gko::make_const_dense_view().
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.
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.
|
inline |
Create a submatrix from the original matrix.
rows | row span |
columns | column span |
size | size of the submatrix (only used for consistency with distributed::Vector) |
|
inline |
Create a submatrix from the original matrix.
rows | row span |
columns | column span |
|
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.
rows | row span |
columns | column span |
stride | stride of the new submatrix. |
Referenced by gko::matrix::Dense< value_type >::create_submatrix().
|
inlinestatic |
Creates a Dense matrix, where the underlying array is a view of another Dense matrix' array.
other | The other matrix on which to create the view |
Referenced by gko::make_dense_view().
|
inlinestatic |
|
inlinestatic |
local_size | Unused |
stride | The stride of the new matrix. |
|
inlinestatic |
stride | The stride of the new matrix. |
|
inlinestatic |
Creates a Dense matrix with the same type as another Dense matrix but on a different executor and with a different size.
other | The other matrix whose type we target. |
exec | The executor of the new matrix. |
size | The size of the new matrix. |
stride | The stride of the new matrix. |
|
overridevirtual |
Extracts the diagonal entries of the matrix into a vector.
diag | the vector into which the diagonal will be written |
Implements gko::DiagonalExtractable< ValueType >.
void gko::matrix::Dense< ValueType >::extract_diagonal | ( | ptr_param< Diagonal< ValueType >> | output | ) | const |
Writes the diagonal of this matrix into an existing diagonal matrix.
output | The output matrix. Its size must match the size of this matrix's diagonal. |
void gko::matrix::Dense< ValueType >::fill | ( | const ValueType | value | ) |
Fill the dense matrix with a given value.
value | the value to be filled |
|
inlinenoexcept |
Returns a pointer to the array of values of the matrix.
|
inlinenoexcept |
Returns the number of elements explicitly stored in the matrix.
|
inlinenoexcept |
Returns the stride of the matrix.
Referenced by gko::matrix::Dense< value_type >::create_submatrix().
|
inlinenoexcept |
Returns a pointer to the array of values of the matrix.
void gko::matrix::Dense< ValueType >::inv_scale | ( | ptr_param< const LinOp > | alpha | ) |
Scales the matrix with the inverse of a scalar.
|
overridevirtual |
Returns a LinOp representing the row permutation of the inverse permuted object.
In the resulting LinOp, the column perm[i]
contains the input column i
.
From the linear algebra perspective, with , this represents the operation
.
permutation_indices | the array of indices containing the permutation order perm . |
Implements gko::Permutable< int32 >.
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.
permutation_indices | The array containing permutation indices. It must have this->get_size()[1] elements. |
output | The output matrix. It must have the dimensions this->get_size() |
|
overridevirtual |
Returns a LinOp representing the row permutation of the inverse permuted object.
In the resulting LinOp, the column perm[i]
contains the input column i
.
From the linear algebra perspective, with , this represents the operation
.
permutation_indices | the array of indices containing the permutation order perm . |
Implements gko::Permutable< int64 >.
void gko::matrix::Dense< ValueType >::inverse_column_permute | ( | const array< int64 > * | permutation_indices, |
ptr_param< Dense< ValueType > > | output | ||
) | 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 , this represents the operation
.
permutation_indices | the array of indices containing the permutation order. |
Reimplemented from gko::Permutable< int32 >.
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.
permutation_indices | The array containing permutation indices. It must have this->get_size()[0] elements. |
output | The output matrix. It must have the dimensions this->get_size() |
|
overridevirtual |
Returns a LinOp representing the symmetric inverse row and column permutation of the Permutable object.
In the resulting LinOp, the entry at location (perm[i],perm[j])
contains the input value (i,j)
.
From the linear algebra perspective, with , this represents the operation
.
permutation_indices | the array of indices containing the permutation order. |
Reimplemented from gko::Permutable< int64 >.
void gko::matrix::Dense< ValueType >::inverse_permute | ( | const array< int64 > * | permutation_indices, |
ptr_param< Dense< ValueType > > | output | ||
) | 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 , this represents the operation
.
permutation_indices | the array of indices containing the permutation order perm . |
Implements gko::Permutable< int32 >.
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.
permutation_indices | The array containing permutation indices. It must have this->get_size()[0] elements. |
output | The output matrix. It must have the dimensions this->get_size() |
|
overridevirtual |
Returns a LinOp representing the row permutation of the inverse permuted object.
In the resulting LinOp, the row perm[i]
contains the input row i
.
From the linear algebra perspective, with , this represents the operation
.
permutation_indices | the array of indices containing the permutation order perm . |
Implements gko::Permutable< int64 >.
void gko::matrix::Dense< ValueType >::inverse_row_permute | ( | const array< int64 > * | permutation_indices, |
ptr_param< Dense< ValueType > > | output | ||
) | const |
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.
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.
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.
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).
|
overridevirtual |
Returns a LinOp representing the symmetric row and column permutation of the Permutable object.
In the resulting LinOp, the entry at location (i,j)
contains the input value (perm[i],perm[j])
.
From the linear algebra perspective, with , this represents the operation
.
permutation_indices | the array of indices containing the permutation order. |
Reimplemented from gko::Permutable< int32 >.
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.
permutation_indices | The array containing permutation indices. It must have this->get_size()[0] elements. |
output | The output matrix. It must have the dimensions this->get_size() |
|
overridevirtual |
Returns a LinOp representing the symmetric row and column permutation of the Permutable object.
In the resulting LinOp, the entry at location (i,j)
contains the input value (perm[i],perm[j])
.
From the linear algebra perspective, with , this represents the operation
.
permutation_indices | the array of indices containing the permutation order. |
Reimplemented from gko::Permutable< int64 >.
void gko::matrix::Dense< ValueType >::permute | ( | const array< int64 > * | permutation_indices, |
ptr_param< Dense< ValueType > > | output | ||
) | const |
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 of this matrix
with the given permutation
.
By default, this computes a symmetric permutation (permute_mode::symmetric). For the effect of the different permutation modes, see permute_mode.
permutation | The input permutation. |
mode | The permutation mode. If permute_mode::inverse is set, we use the inverse permutation ![]() ![]() |
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.
output | the output matrix. |
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 of this matrix
with the given row and column permutations
and
.
The operation will compute , or
if
invert
is false
, and , or
if
invert
is true
.
row_permutation | The permutation ![]() |
column_permutation | The permutation ![]() |
invert | If set to false , uses the input permutations, otherwise uses their inverses ![]() |
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 |
std::unique_ptr<Dense> gko::matrix::Dense< ValueType >::permute | ( | ptr_param< const Permutation< int64 >> | permutation, |
permute_mode | mode = permute_mode::symmetric |
||
) | const |
void gko::matrix::Dense< ValueType >::permute | ( | ptr_param< const Permutation< int64 >> | permutation, |
ptr_param< Dense< ValueType > > | output, | ||
permute_mode | mode | ||
) | const |
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 |
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 |
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.
gather_indices | pointer to an array containing row indices from this matrix. It may contain duplicates. |
gather_indices->get_size()
rows containing the gathered rows from this matrix: output(i,j) = input(gather_indices(i), j)
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
gather_indices | pointer to an array containing row indices from this matrix. It may contain duplicates. |
row_collection | pointer 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_size() rows. |
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.
gather_indices | pointer to an array containing row indices from this matrix. It may contain duplicates. |
gather_indices->get_size()
rows containing the gathered rows from this matrix: output(i,j) = input(gather_indices(i), j)
void gko::matrix::Dense< ValueType >::row_gather | ( | const array< int64 > * | gather_indices, |
ptr_param< LinOp > | row_collection | ||
) | const |
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.
alpha | scaling the result of row gathering |
gather_indices | pointer to an array containing row indices from this matrix. It may contain duplicates. |
beta | scaling the input row_collection |
row_collection | pointer 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_size() rows. |
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 |
|
overridevirtual |
Returns a LinOp representing the row permutation of the Permutable object.
In the resulting LinOp, the row i
contains the input row perm[i]
.
From the linear algebra perspective, with , this represents the operation
.
permutation_indices | the array of indices containing the permutation order. |
Implements gko::Permutable< int32 >.
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.
permutation_indices | The array containing permutation indices. It must have this->get_size()[0] elements. |
output | The output matrix. It must have the dimensions this->get_size() |
|
overridevirtual |
Returns a LinOp representing the row permutation of the Permutable object.
In the resulting LinOp, the row i
contains the input row perm[i]
.
From the linear algebra perspective, with , this represents the operation
.
permutation_indices | the array of indices containing the permutation order. |
Implements gko::Permutable< int64 >.
void gko::matrix::Dense< ValueType >::row_permute | ( | const array< int64 > * | permutation_indices, |
ptr_param< Dense< ValueType > > | output | ||
) | const |
void gko::matrix::Dense< ValueType >::scale | ( | ptr_param< const LinOp > | alpha | ) |
Scales the matrix with a scalar (aka: BLAS scal).
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)
permutation | The scaled permutation. |
mode | The permutation mode. |
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.
output | the output matrix. |
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)
row_permutation | The scaled row permutation. |
column_permutation | The scaled column permutation. |
invert | If set to false , uses the input permutations, otherwise uses their inverses ![]() |
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.
output | the output matrix. |
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 |
void gko::matrix::Dense< ValueType >::scale_permute | ( | ptr_param< const ScaledPermutation< value_type, int64 >> | permutation, |
ptr_param< Dense< ValueType > > | output, | ||
permute_mode | mode | ||
) | const |
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 |
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 |
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).
|
overridevirtual |
Returns a LinOp representing the transpose of the Transposable object.
Implements gko::Transposable.
void gko::matrix::Dense< ValueType >::transpose | ( | ptr_param< Dense< ValueType > > | output | ) | const |
Writes the transposed matrix into the given output matrix.
output | The output matrix. It must have the dimensions gko::transpose(this->get_size()) |