Morpheus 1.0.0
Dynamic matrix type and algorithms for sparse matrices
Loading...
Searching...
No Matches
Public Types | Public Member Functions | List of all members
Morpheus::DiaMatrix< ValueType, Properties > Class Template Reference

Implementation of the Diagonal (DIA) Sparse Matrix Format Representation. More...

#include <Morpheus_DiaMatrix.hpp>

Inheritance diagram for Morpheus::DiaMatrix< ValueType, Properties >:
Inheritance graph
[legend]
Collaboration diagram for Morpheus::DiaMatrix< ValueType, Properties >:
Collaboration graph
[legend]

Public Types

using traits = ContainerTraits< DiaMatrix, ValueType, Properties... >
 
using type = typename traits::type
 
using base = MatrixBase< DiaMatrix, ValueType, Properties... >
 
using tag = typename MatrixFormatTag< DiaFormatTag >::tag
 
using value_type = typename traits::value_type
 
using non_const_value_type = typename traits::non_const_value_type
 
using size_type = typename traits::size_type
 
using index_type = typename traits::index_type
 
using non_const_index_type = typename traits::non_const_index_type
 
using array_layout = typename traits::array_layout
 
using backend = typename traits::backend
 
using memory_space = typename traits::memory_space
 
using execution_space = typename traits::execution_space
 
using device_type = typename traits::device_type
 
using memory_traits = typename traits::memory_traits
 
using HostMirror = typename traits::HostMirror
 
using pointer = typename traits::pointer
 
using const_pointer = typename traits::const_pointer
 
using reference = typename traits::reference
 
using const_reference = typename traits::const_reference
 
using index_array_type = Morpheus::DenseVector< index_type, size_type, array_layout, backend, memory_traits >
 
using const_index_array_type = const index_array_type
 
using index_array_pointer = typename index_array_type::value_array_pointer
 
using index_array_reference = typename index_array_type::value_array_reference
 
using const_index_array_reference = const index_array_reference
 
using value_array_type = Morpheus::DenseMatrix< value_type, size_type, array_layout, backend, memory_traits >
 
using const_value_array_type = const value_array_type
 
using value_array_pointer = typename value_array_type::value_array_pointer
 
using value_array_reference = typename value_array_type::value_array_reference
 
using const_value_array_reference = const value_array_reference
 
- Public Types inherited from Morpheus::MatrixBase< DiaMatrix, ValueType, Properties... >
using type = MatrixBase< DiaMatrix, ValueType, Properties... >
 The traits associated with the particular container.
 
using traits = ContainerTraits< DiaMatrix, ValueType, Properties... >
 The type of the indices held by the container.
 
using size_type = typename traits::size_type
 
- Public Types inherited from Morpheus::ContainerTraits< Container, ValueType, Properties >
enum  { is_hostspace = std::is_same<MemorySpace, Kokkos::HostSpace>::value }
 
enum  { is_managed = MemoryTraits::is_unmanaged == 0 }
 
using value_type = ValueType
 The type of values held by the container.
 
using const_value_type = typename std::add_const< ValueType >::type
 The const type of values held by the container.
 
using non_const_value_type = typename std::remove_const< ValueType >::type
 The non-const type of values held by the container.
 
using index_type = IndexType
 The type of indices held by the container.
 
using size_type = size_t
 The size type of the container.
 
using non_const_index_type = typename std::remove_const< IndexType >::type
 The non-const type of indices held by the container.
 
using array_layout = ArrayLayout
 The storage layout of data held by the container.
 
using backend = Backend
 The backend out of which algorithms will be dispatched from.
 
using execution_space = ExecutionSpace
 The space in which member functions will be executed in.
 
using memory_space = MemorySpace
 The space in which data will be stored in.
 
using device_type = Morpheus::Device< execution_space, memory_space, backend >
 A device aware of the execution, memory spaces and backend.
 
using memory_traits = MemoryTraits
 Represents the user's intended access behaviour.
 
using host_mirror_backend = typename Morpheus::HostMirror< backend >::backend
 The host equivalent backend.
 
using type = Container< value_type, index_type, array_layout, backend, memory_traits >
 The complete type of the container.
 
using HostMirror = Container< non_const_value_type, non_const_index_type, array_layout, Morpheus::Device< typename host_mirror_backend::execution_space, typename host_mirror_backend::memory_space, typename host_mirror_backend::backend >, typename Kokkos::MemoryManaged >
 The host mirror equivalent for the container. More...
 
using pointer = typename std::add_pointer< type >::type
 The pointer type of the container.
 
using const_pointer = typename std::add_pointer< typename std::add_const< type >::type >::type
 The const pointer type of the container.
 
using reference = typename std::add_lvalue_reference< type >::type
 The reference type of the container.
 
using const_reference = typename std::add_lvalue_reference< typename std::add_const< type >::type >::type
 The const reference type of the container.
 

Public Member Functions

 ~DiaMatrix ()=default
 The default destructor.
 
 DiaMatrix (const DiaMatrix &)=default
 The default copy contructor (shallow copy) of a DiaMatrix container from another DiaMatrix container with the same properties.
 
 DiaMatrix (DiaMatrix &&)=default
 The default move contructor (shallow copy) of a DiaMatrix container from another DiaMatrix container with the same properties.
 
DiaMatrixoperator= (const DiaMatrix &)=default
 The default copy assignment (shallow copy) of a DiaMatrix container from another DiaMatrix container with the same properties.
 
DiaMatrixoperator= (DiaMatrix &&)=default
 The default move assignment (shallow copy) of a DiaMatrix container from another DiaMatrix container with the same properties.
 
 DiaMatrix ()
 Construct an empty DiaMatrix object.
 
 DiaMatrix (const size_type num_rows, const size_type num_cols, const size_type num_entries, const size_type num_diagonals, const size_type alignment=32)
 Construct a DiaMatrix object with shape (num_rows, num_cols) and number of non-zeros equal to num_entries concentrated around num_diagonals diagonals. More...
 
template<typename ValueArray , typename IndexArray >
 DiaMatrix (const size_type num_rows, const size_type num_cols, const size_type num_entries, const IndexArray &diag_offsets, const ValueArray &vals, typename std::enable_if< is_dense_matrix_format_container< ValueArray >::value &&is_dense_vector_format_container< IndexArray >::value &&is_compatible< typename DiaMatrix::value_array_type, ValueArray >::value &&is_compatible< typename DiaMatrix::index_array_type, IndexArray >::value &&!ValueArray::memory_traits::is_unmanaged &&!IndexArray::memory_traits::is_unmanaged >::type *=nullptr)
 Construct a DiaMatrix object with shape (num_rows, num_cols) and number of non-zeros equal to num_entries and assign the diagonal offsets and values from DenseVector and DenseMatrix arrays equivalently. More...
 
template<typename ValuePtr , typename IndexPtr >
 DiaMatrix (const size_type num_rows, const size_type num_cols, const size_type num_entries, IndexPtr diag_offsets_ptr, ValuePtr vals_ptr, const size_type num_diagonals, const size_type alignment=32, typename std::enable_if<(std::is_pointer< ValuePtr >::value &&is_same_value_type< value_type, ValuePtr >::value &&memory_traits::is_unmanaged) &&(std::is_pointer< IndexPtr >::value &&is_same_index_type< index_type, IndexPtr >::value &&memory_traits::is_unmanaged)>::type *=nullptr)
 
template<class VR , class... PR>
 DiaMatrix (const DiaMatrix< VR, PR... > &src, typename std::enable_if< is_format_compatible< DiaMatrix, DiaMatrix< VR, PR... > >::value >::type *=nullptr)
 Constructs a DiaMatrix from another compatible DiaMatrix. More...
 
template<class VR , class... PR>
std::enable_if< is_format_compatible< DiaMatrix, DiaMatrix< VR, PR... > >::value, DiaMatrix & >::type operator= (const DiaMatrix< VR, PR... > &src)
 Assigns a DiaMatrix from another compatible DiaMatrix. More...
 
template<class VR , class... PR>
 DiaMatrix (const DynamicMatrix< VR, PR... > &src, typename std::enable_if< is_dynamically_compatible< DiaMatrix, DynamicMatrix< VR, PR... > >::value >::type *=nullptr)
 Constructs a DiaMatrix from a compatible DynamicMatrix. More...
 
template<class VR , class... PR>
std::enable_if< is_dynamically_compatible< DiaMatrix, DynamicMatrix< VR, PR... > >::value, DiaMatrix & >::type operator= (const DynamicMatrix< VR, PR... > &src)
 Assigns a DiarMatrix from a compatible DynamicMatrix. More...
 
template<typename MatrixType >
 DiaMatrix (const MatrixType &src)=delete
 Construct a DiaMatrix object from another storage format. This functionality is disabled to avoid implicit copies and conversion operations. More...
 
template<typename MatrixType >
reference operator= (const MatrixType &src)=delete
 Assigns to DiaMatrix object from another storage format. This functionality is disabled to avoid implicit copies and conversion operations. More...
 
void resize (const size_type num_rows, const size_type num_cols, const size_type num_entries, const size_type num_diagonals, const size_type alignment=32)
 Resizes DiaMatrix with shape of (num_rows, num_cols) and sets number of non-zero entries to num_entries, number of diagonals to num_diagonals and the alignment. More...
 
template<class VR , class... PR>
void resize (const DiaMatrix< VR, PR... > &src)
 Resizes DiaMatrix with the shape and number of non-zero entries of another DiaMatrix with different parameters. More...
 
template<class VR , class... PR>
DiaMatrixallocate (const DiaMatrix< VR, PR... > &src)
 Allocates memory from another DiaMatrix container with different properties. More...
 
formats_e format_enum () const
 Returns the format enum assigned to the DiaMatrix container. More...
 
int format_index () const
 Returns the equivalent index to the format enum assigned to the DiaMatrix container. More...
 
MORPHEUS_FORCEINLINE_FUNCTION index_array_reference diagonal_offsets (size_type n)
 Returns a reference to the diagonal offsets of the matrix with index n. More...
 
MORPHEUS_FORCEINLINE_FUNCTION value_array_reference values (size_type i, size_type j)
 Returns a reference to the value of the matrix with indexes (i, j) More...
 
MORPHEUS_FORCEINLINE_FUNCTION const_index_array_reference cdiagonal_offsets (size_type n) const
 Returns a const-reference to the diagonal offsets of the matrix with index n. More...
 
MORPHEUS_FORCEINLINE_FUNCTION const_value_array_reference cvalues (size_type i, size_type j) const
 Returns a const-reference to the value of the matrix with indexes (i, j) More...
 
MORPHEUS_FORCEINLINE_FUNCTION index_array_typediagonal_offsets ()
 Returns a reference to the diagonal offsets of the matrix. More...
 
MORPHEUS_FORCEINLINE_FUNCTION value_array_typevalues ()
 Returns a reference to the values of the matrix. More...
 
MORPHEUS_FORCEINLINE_FUNCTION const_index_array_typecdiagonal_offsets () const
 Returns a const-reference to the diagonal offsets of the matrix. More...
 
MORPHEUS_FORCEINLINE_FUNCTION const_value_array_typecvalues () const
 Returns a const-reference to the values of the matrix. More...
 
MORPHEUS_FORCEINLINE_FUNCTION size_type ndiags () const
 Returns the number of occupied diagonals of the matrix. More...
 
MORPHEUS_FORCEINLINE_FUNCTION size_type alignment () const
 Returns the amount of padding used to align the data. More...
 
MORPHEUS_FORCEINLINE_FUNCTION void set_ndiags (const size_type num_diagonals)
 Sets the number of occupied diagonals of the matrix. More...
 
MORPHEUS_FORCEINLINE_FUNCTION void set_alignment (const size_type alignment)
 Sets amount of padding with which to align the data. More...
 
- Public Member Functions inherited from Morpheus::MatrixBase< DiaMatrix, ValueType, Properties... >
 MatrixBase ()
 Default constructor. More...
 
 MatrixBase (size_type rows, size_type cols, size_type entries=0)
 Construct a MatrixBase object with shape (num_rows, num_cols) and number of non-zeros equal to num_entries. More...
 
void resize (size_type rows, size_type cols, size_type entries)
 Resizes MatrixBase with shape of (num_rows, num_cols) and sets number of non-zero entries to num_entries. More...
 
size_type nrows () const
 Number of rows of the matrix. More...
 
size_type ncols () const
 Number of columns of the matrix. More...
 
size_type nnnz () const
 Number of non-zeros of the matrix. More...
 
void set_nrows (const size_type rows)
 Set the number of rows of the matrix. More...
 
void set_ncols (const size_type cols)
 Set the number of columns of the matrix. More...
 
void set_nnnz (const size_type nnz)
 Set the number of non-zeros of the matrix. More...
 
MatrixStructure structure () const
 The specialized structure of the matrix e.g Symmetric. More...
 
MatrixOptions options () const
 Information about specific characteristics of the matrix e.g has short rows. More...
 
void set_structure (MatrixStructure op)
 Set the structure of the matrix. More...
 
void set_options (MatrixOptions op)
 Set the characteristics of the matrix. More...
 

Detailed Description

template<class ValueType, class... Properties>
class Morpheus::DiaMatrix< ValueType, Properties >

Implementation of the Diagonal (DIA) Sparse Matrix Format Representation.

Template Parameters
ValueTypeType of values to store
PropertiesOptional properties to modify the behaviour of the container. Sensible defaults are selected based on the configuration. Please refer to impl/Morpheus_ContainerTraits.hpp to find out more about the valid properties.
Overview
The DiaMatrix container is a two-dimensional container that represents a sparse matrix. This container stores the contents of a sparse matrix along its diagonals allowing for very quick arithmetic operations for matrices with very few non-zero diagonals. However, in the case if the number of non-zero diagonals becomes very large the performance is hindered because of the excessive zero padding required. It is a polymorphic container in the sense that it can store scalar or integer type values, on host or device depending how the template parameters are selected.
Example
#include <Morpheus_Core.hpp>
// Matrix to Build
// [1 * 2]
// [* * 3]
// [* 4 *]
int main(){
using index_array_type = typename Matrix::index_array_type;
using value_array_type = typename Matrix::value_array_type;
index_array_type offsets(4, 0);
offset[0] = -1; offset[1] = 0; offset[2] = 1; offset[3] = 2;
// [* 1 0 2]
// [0 0 3 *]
// [4 0 * *]
value(0,0) = 0; value(1,0) = 0; value(2,0) = 4;
value(0,1) = 1; value(1,1) = 0; value(2,1) = 0;
value(0,2) = 0; value(1,2) = 3; value(2,2) = 0;
value(0,3) = 2; value(1,3) = 0; value(2,3) = 0;
// Construct the matrix from diagonal offsets and values
Matrix A(3, 3, 4, offset, values);
Morpheus::print(A); // prints A
}
Implementation of the Diagonal (DIA) Sparse Matrix Format Representation.
Definition: Morpheus_DiaMatrix.hpp:97
MORPHEUS_FORCEINLINE_FUNCTION value_array_type & values()
Returns a reference to the values of the matrix.
Definition: Morpheus_DiaMatrix.hpp:510

Member Typedef Documentation

◆ index_array_type

template<class ValueType , class... Properties>
using Morpheus::DiaMatrix< ValueType, Properties >::index_array_type = Morpheus::DenseVector<index_type, size_type, array_layout, backend, memory_traits>

The type of DenseVector that holds the index_type data

◆ index_type

template<class ValueType , class... Properties>
using Morpheus::DiaMatrix< ValueType, Properties >::index_type = typename traits::index_type

The type of the indices held by the container - can be const

◆ non_const_index_type

template<class ValueType , class... Properties>
using Morpheus::DiaMatrix< ValueType, Properties >::non_const_index_type = typename traits::non_const_index_type

The non-constant type of the indices held by the container

◆ non_const_value_type

template<class ValueType , class... Properties>
using Morpheus::DiaMatrix< ValueType, Properties >::non_const_value_type = typename traits::non_const_value_type

The non-constant type of the values held by the container

◆ tag

template<class ValueType , class... Properties>
using Morpheus::DiaMatrix< ValueType, Properties >::tag = typename MatrixFormatTag<DiaFormatTag>::tag

The tag associated specificaly to the particular container

◆ traits

template<class ValueType , class... Properties>
using Morpheus::DiaMatrix< ValueType, Properties >::traits = ContainerTraits<DiaMatrix, ValueType, Properties...>

The traits associated with the particular container

◆ type

template<class ValueType , class... Properties>
using Morpheus::DiaMatrix< ValueType, Properties >::type = typename traits::type

The complete type of the container

◆ value_array_type

template<class ValueType , class... Properties>
using Morpheus::DiaMatrix< ValueType, Properties >::value_array_type = Morpheus::DenseMatrix<value_type, size_type, array_layout, backend, memory_traits>

The type of DenseMatrix that holds the value_type data

◆ value_type

template<class ValueType , class... Properties>
using Morpheus::DiaMatrix< ValueType, Properties >::value_type = typename traits::value_type

The type of the values held by the container - can be const

Constructor & Destructor Documentation

◆ DiaMatrix() [1/5]

template<class ValueType , class... Properties>
Morpheus::DiaMatrix< ValueType, Properties >::DiaMatrix ( const size_type  num_rows,
const size_type  num_cols,
const size_type  num_entries,
const size_type  num_diagonals,
const size_type  alignment = 32 
)
inline

Construct a DiaMatrix object with shape (num_rows, num_cols) and number of non-zeros equal to num_entries concentrated around num_diagonals diagonals.

Parameters
num_rowsNumber of rows of the matrix.
num_colsNumber of columns of the matrix.
num_entriesNumber of non-zero values in the matrix.
num_diagonalsNumber of occupied diagonals
alignmentAmount of padding used to align the data (default=32)

◆ DiaMatrix() [2/5]

template<class ValueType , class... Properties>
template<typename ValueArray , typename IndexArray >
Morpheus::DiaMatrix< ValueType, Properties >::DiaMatrix ( const size_type  num_rows,
const size_type  num_cols,
const size_type  num_entries,
const IndexArray &  diag_offsets,
const ValueArray &  vals,
typename std::enable_if< is_dense_matrix_format_container< ValueArray >::value &&is_dense_vector_format_container< IndexArray >::value &&is_compatible< typename DiaMatrix< ValueType, Properties >::value_array_type, ValueArray >::value &&is_compatible< typename DiaMatrix< ValueType, Properties >::index_array_type, IndexArray >::value &&!ValueArray::memory_traits::is_unmanaged &&!IndexArray::memory_traits::is_unmanaged >::type = nullptr 
)
inlineexplicit

Construct a DiaMatrix object with shape (num_rows, num_cols) and number of non-zeros equal to num_entries and assign the diagonal offsets and values from DenseVector and DenseMatrix arrays equivalently.

Template Parameters
ValueArrayValue type DenseMatrix type.
IndexArrayIndex type DenseVector type.
Parameters
num_rowsNumber of rows of the matrix.
num_colsNumber of columns of the matrix.
num_entriesNumber of non-zero values in the matrix.
diag_offsetsDenseVector containing the diagonal offsets of the matrix.
valsDenseMatrix containing the values of the matrix.

◆ DiaMatrix() [3/5]

template<class ValueType , class... Properties>
template<class VR , class... PR>
Morpheus::DiaMatrix< ValueType, Properties >::DiaMatrix ( const DiaMatrix< VR, PR... > &  src,
typename std::enable_if< is_format_compatible< DiaMatrix< ValueType, Properties >, DiaMatrix< VR, PR... > >::value >::type = nullptr 
)
inline

Constructs a DiaMatrix from another compatible DiaMatrix.

Constructs a DiaMatrix from another compatible DiaMatrix i.e a
matrix that satisfies the is_format_compatible check.
Template Parameters
VRType of Values the Other Matrix holds.
PRProperties of the Other Matrix.
Parameters
srcThe matrix we are constructing from.

◆ DiaMatrix() [4/5]

template<class ValueType , class... Properties>
template<class VR , class... PR>
Morpheus::DiaMatrix< ValueType, Properties >::DiaMatrix ( const DynamicMatrix< VR, PR... > &  src,
typename std::enable_if< is_dynamically_compatible< DiaMatrix< ValueType, Properties >, DynamicMatrix< VR, PR... > >::value >::type = nullptr 
)
inline

Constructs a DiaMatrix from a compatible DynamicMatrix.

Overview
Constructs a DiaMatrix from a compatible DynamicMatrix i.e a matrix that satisfies the is_dynamically_compatible check. Note that when the active type of the dynamic matrix is different from the concrete type, this will result in an exception thrown.
Template Parameters
VRType of Values the Other Matrix holds.
PRProperties of the Other Matrix.
Parameters
srcThe matrix we are constructing from.

◆ DiaMatrix() [5/5]

template<class ValueType , class... Properties>
template<typename MatrixType >
Morpheus::DiaMatrix< ValueType, Properties >::DiaMatrix ( const MatrixType &  src)
delete

Construct a DiaMatrix object from another storage format. This functionality is disabled to avoid implicit copies and conversion operations.

Template Parameters
MatrixTypeAny of the supported storage formats.
Parameters
srcThe source container.

Member Function Documentation

◆ alignment()

template<class ValueType , class... Properties>
MORPHEUS_FORCEINLINE_FUNCTION size_type Morpheus::DiaMatrix< ValueType, Properties >::alignment ( ) const
inline

Returns the amount of padding used to align the data.

Returns
size_type Amount of padding used to align the data

◆ allocate()

template<class ValueType , class... Properties>
template<class VR , class... PR>
DiaMatrix & Morpheus::DiaMatrix< ValueType, Properties >::allocate ( const DiaMatrix< VR, PR... > &  src)
inline

Allocates memory from another DiaMatrix container with different properties.

Template Parameters
VRValue Type of the container we are allocating from.
PROptional properties of the container we are allocating from.
Parameters
srcThe DiaMatrix container we are allocating from.

◆ cdiagonal_offsets() [1/2]

template<class ValueType , class... Properties>
MORPHEUS_FORCEINLINE_FUNCTION const_index_array_type & Morpheus::DiaMatrix< ValueType, Properties >::cdiagonal_offsets ( ) const
inline

Returns a const-reference to the diagonal offsets of the matrix.

Returns
index_array_type& A reference to the diagonal offsets.

◆ cdiagonal_offsets() [2/2]

template<class ValueType , class... Properties>
MORPHEUS_FORCEINLINE_FUNCTION const_index_array_reference Morpheus::DiaMatrix< ValueType, Properties >::cdiagonal_offsets ( size_type  n) const
inline

Returns a const-reference to the diagonal offsets of the matrix with index n.

Parameters
nIndex of the value to extract
Returns
Diagonal offset at index n

◆ cvalues() [1/2]

template<class ValueType , class... Properties>
MORPHEUS_FORCEINLINE_FUNCTION const_value_array_type & Morpheus::DiaMatrix< ValueType, Properties >::cvalues ( ) const
inline

Returns a const-reference to the values of the matrix.

Returns
value_array_type& A reference to the values.

◆ cvalues() [2/2]

template<class ValueType , class... Properties>
MORPHEUS_FORCEINLINE_FUNCTION const_value_array_reference Morpheus::DiaMatrix< ValueType, Properties >::cvalues ( size_type  i,
size_type  j 
) const
inline

Returns a const-reference to the value of the matrix with indexes (i, j)

Parameters
iRow index of the value to extract
jColumn index of the value to extract
Returns
Value of the element at index n

◆ diagonal_offsets() [1/2]

template<class ValueType , class... Properties>
MORPHEUS_FORCEINLINE_FUNCTION index_array_type & Morpheus::DiaMatrix< ValueType, Properties >::diagonal_offsets ( )
inline

Returns a reference to the diagonal offsets of the matrix.

Returns
index_array_type& A reference to the diagonal offsets.

◆ diagonal_offsets() [2/2]

template<class ValueType , class... Properties>
MORPHEUS_FORCEINLINE_FUNCTION index_array_reference Morpheus::DiaMatrix< ValueType, Properties >::diagonal_offsets ( size_type  n)
inline

Returns a reference to the diagonal offsets of the matrix with index n.

Parameters
nIndex of the value to extract
Returns
Diagonal offset at index n

◆ format_enum()

template<class ValueType , class... Properties>
formats_e Morpheus::DiaMatrix< ValueType, Properties >::format_enum ( ) const
inline

Returns the format enum assigned to the DiaMatrix container.

Returns
formats_e The format enum

◆ format_index()

template<class ValueType , class... Properties>
int Morpheus::DiaMatrix< ValueType, Properties >::format_index ( ) const
inline

Returns the equivalent index to the format enum assigned to the DiaMatrix container.

Returns
int The equivalent index to format_e

◆ ndiags()

template<class ValueType , class... Properties>
MORPHEUS_FORCEINLINE_FUNCTION size_type Morpheus::DiaMatrix< ValueType, Properties >::ndiags ( ) const
inline

Returns the number of occupied diagonals of the matrix.

Returns
size_type Number of occupied diagonals

◆ operator=() [1/3]

template<class ValueType , class... Properties>
template<class VR , class... PR>
std::enable_if< is_format_compatible< DiaMatrix, DiaMatrix< VR, PR... > >::value, DiaMatrix & >::type Morpheus::DiaMatrix< ValueType, Properties >::operator= ( const DiaMatrix< VR, PR... > &  src)
inline

Assigns a DiaMatrix from another compatible DiaMatrix.

Overview
Assigns a DiaMatrix from another compatible DiaMatrix i.e a matrix that satisfies the is_format_compatible check.
Template Parameters
VRType of Values the Other Matrix holds.
PRProperties of the Other Matrix.
Parameters
srcThe matrix we are assigning from.

◆ operator=() [2/3]

template<class ValueType , class... Properties>
template<class VR , class... PR>
std::enable_if< is_dynamically_compatible< DiaMatrix, DynamicMatrix< VR, PR... > >::value, DiaMatrix & >::type Morpheus::DiaMatrix< ValueType, Properties >::operator= ( const DynamicMatrix< VR, PR... > &  src)
inline

Assigns a DiarMatrix from a compatible DynamicMatrix.

Overview
Assigns a DiarMatrix from a compatible DynamicMatrix i.e a matrix that satisfies the is_dynamically_compatible check. Note that when the active type of the dynamic matrix is different from the concrete type, this will result in an exception thrown.
Template Parameters
VRType of Values the Other Matrix holds.
PRProperties of the Other Matrix.
Parameters
srcThe matrix we are assigning from.

◆ operator=() [3/3]

template<class ValueType , class... Properties>
template<typename MatrixType >
reference Morpheus::DiaMatrix< ValueType, Properties >::operator= ( const MatrixType &  src)
delete

Assigns to DiaMatrix object from another storage format. This functionality is disabled to avoid implicit copies and conversion operations.

Template Parameters
MatrixTypeAny of the supported storage formats.
Parameters
srcThe source container.

◆ resize() [1/2]

template<class ValueType , class... Properties>
template<class VR , class... PR>
void Morpheus::DiaMatrix< ValueType, Properties >::resize ( const DiaMatrix< VR, PR... > &  src)
inline

Resizes DiaMatrix with the shape and number of non-zero entries of another DiaMatrix with different parameters.

Template Parameters
VRType of values the source matrix stores.
PROther properties of source matrix.
Parameters
srcThe source CsrMatrix we are resizing from.

◆ resize() [2/2]

template<class ValueType , class... Properties>
void Morpheus::DiaMatrix< ValueType, Properties >::resize ( const size_type  num_rows,
const size_type  num_cols,
const size_type  num_entries,
const size_type  num_diagonals,
const size_type  alignment = 32 
)
inline

Resizes DiaMatrix with shape of (num_rows, num_cols) and sets number of non-zero entries to num_entries, number of diagonals to num_diagonals and the alignment.

Parameters
num_rowsNumber of rows of resized matrix.
num_colsNumber of columns of resized matrix.
num_entriesNumber of non-zero entries in resized matrix.
num_diagonalsNumber of occupied diagonals in resized matrix.
alignmentAmount of padding used to align the data (default=32)

◆ set_alignment()

template<class ValueType , class... Properties>
MORPHEUS_FORCEINLINE_FUNCTION void Morpheus::DiaMatrix< ValueType, Properties >::set_alignment ( const size_type  alignment)
inline

Sets amount of padding with which to align the data.

Parameters
alignmentNew amount of padding.

◆ set_ndiags()

template<class ValueType , class... Properties>
MORPHEUS_FORCEINLINE_FUNCTION void Morpheus::DiaMatrix< ValueType, Properties >::set_ndiags ( const size_type  num_diagonals)
inline

Sets the number of occupied diagonals of the matrix.

Parameters
num_diagonalsnumber of new diagonals

◆ values() [1/2]

template<class ValueType , class... Properties>
MORPHEUS_FORCEINLINE_FUNCTION value_array_type & Morpheus::DiaMatrix< ValueType, Properties >::values ( )
inline

Returns a reference to the values of the matrix.

Returns
value_array_type& A reference to the values.

◆ values() [2/2]

template<class ValueType , class... Properties>
MORPHEUS_FORCEINLINE_FUNCTION value_array_reference Morpheus::DiaMatrix< ValueType, Properties >::values ( size_type  i,
size_type  j 
)
inline

Returns a reference to the value of the matrix with indexes (i, j)

Parameters
iRow index of the value to extract
jColumn index of the value to extract
Returns
Value of the element at index n

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