Template Numerical Library version\ main:94209208
Loading...
Searching...
No Matches
Public Types | Public Member Functions | List of all members
TNL::Matrices::MultidiagonalMatrixView< Real, Device, Index, Organization > Class Template Reference

Implementation of sparse multidiagonal matrix. More...

#include <TNL/Matrices/MultidiagonalMatrixView.h>

Inheritance diagram for TNL::Matrices::MultidiagonalMatrixView< Real, Device, Index, Organization >:
Inheritance graph
[legend]
Collaboration diagram for TNL::Matrices::MultidiagonalMatrixView< Real, Device, Index, Organization >:
Collaboration graph
[legend]

Public Types

using ConstViewType = MultidiagonalMatrixView< std::add_const_t< Real >, Device, Index, Organization >
 Matrix view type for constant instances.
 
using DeviceType = Device
 The device where the matrix is allocated.
 
using IndexType = Index
 The type used for matrix elements indexing.
 
using RealType = typename Base::RealType
 The type of matrix elements.
 
template<typename _Real = Real, typename _Device = Device, typename _Index = Index, ElementsOrganization Organization_ = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization()>
using Self = MultidiagonalMatrixView< _Real, _Device, _Index, Organization_ >
 Helper type for getting self type or its modifications.
 
using ViewType = MultidiagonalMatrixView< Real, Device, Index, Organization >
 Type of related matrix view.
 
- Public Types inherited from TNL::Matrices::MultidiagonalMatrixBase< Real, Device, Index, Organization >
using ConstRowView = typename RowView::ConstRowView
 Type for accessing constant matrix rows.
 
using DeviceType = Device
 The device where the matrix is allocated.
 
using DiagonalOffsetsView
 
using HostDiagonalOffsetsView
 
using IndexerType = details::MultidiagonalMatrixIndexer< Index, Organization == Algorithms::Segments::RowMajorOrder >
 
using IndexType = Index
 The type used for matrix elements indexing.
 
using RealType = typename Base::RealType
 The type of matrix elements.
 
using RowView = MultidiagonalMatrixRowView< typename Base::ValuesViewType, IndexerType, DiagonalOffsetsView >
 Type for accessing matrix rows.
 
- Public Types inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
using ConstValuesViewType
 Type of constant vector view holding values of matrix elements.
 
using DeviceType
 The device where the matrix is allocated.
 
using IndexType
 The type used for matrix elements indexing.
 
using RealType
 The type of matrix elements.
 
using RowCapacitiesType
 
using ValuesViewType
 Type of vector view holding values of matrix elements.
 

Public Member Functions

__cuda_callable__ MultidiagonalMatrixView ()=default
 Constructor with no parameters.
 
__cuda_callable__ MultidiagonalMatrixView (const MultidiagonalMatrixView &)=default
 Copy constructor.
 
__cuda_callable__ MultidiagonalMatrixView (MultidiagonalMatrixView &&) noexcept=default
 Move constructor.
 
__cuda_callable__ MultidiagonalMatrixView (typename Base::ValuesViewType values, typename Base::DiagonalOffsetsView diagonalOffsets, typename Base::HostDiagonalOffsetsView hostDiagonalOffsets, typename Base::IndexerType indexer)
 Constructor with all necessary data and views.
 
__cuda_callable__ void bind (MultidiagonalMatrixView view)
 Method for rebinding (reinitialization) using another multidiagonal matrix view.
 
ConstViewType getConstView () const
 Returns a non-modifiable view of the multidiagonal matrix.
 
ViewType getView ()
 Returns a modifiable view of the multidiagonal matrix.
 
MultidiagonalMatrixViewoperator= (const MultidiagonalMatrixView &)=delete
 Copy-assignment operator.
 
MultidiagonalMatrixViewoperator= (MultidiagonalMatrixView &&)=delete
 Move-assignment operator.
 
void save (File &file) const
 Method for saving the matrix to a file.
 
- Public Member Functions inherited from TNL::Matrices::MultidiagonalMatrixBase< Real, Device, Index, Organization >
__cuda_callable__ MultidiagonalMatrixBase ()=default
 Constructor with no parameters.
 
__cuda_callable__ MultidiagonalMatrixBase (const MultidiagonalMatrixBase &)=default
 Copy constructor.
 
__cuda_callable__ MultidiagonalMatrixBase (MultidiagonalMatrixBase &&) noexcept=default
 Move constructor.
 
__cuda_callable__ MultidiagonalMatrixBase (typename Base::ValuesViewType values, DiagonalOffsetsView diagonalOffsets, HostDiagonalOffsetsView hostDiagonalOffsets, IndexerType indexer)
 Constructor with all necessary data and views.
 
__cuda_callable__ void addElement (IndexType row, IndexType column, const RealType &value, const RealType &thisElementMultiplicator=1.0)
 Add element at given row and column to given value.
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_>
void addMatrix (const MultidiagonalMatrixBase< Real_, Device_, Index_, Organization_ > &matrix, const RealType &matrixMultiplicator=1.0, const RealType &thisMatrixMultiplicator=1.0)
 
template<typename Function >
void forAllElements (Function &function)
 This method calls forElements for all matrix rows.
 
template<typename Function >
void forAllElements (Function &function) const
 This method calls forElements for all matrix rows (for constant instances).
 
template<typename Function >
void forAllRows (Function &&function)
 Method for parallel iteration over all matrix rows.
 
template<typename Function >
void forAllRows (Function &&function) const
 Method for parallel iteration over all matrix rows for constant instances.
 
template<typename Function >
void forElements (IndexType begin, IndexType end, Function &function)
 Method for iteration over all matrix rows for non-constant instances.
 
template<typename Function >
void forElements (IndexType begin, IndexType end, Function &function) const
 Method for iteration over all matrix rows for constant instances.
 
template<typename Function >
void forRows (IndexType begin, IndexType end, Function &&function)
 Method for parallel iteration over matrix rows from interval [begin, end).
 
template<typename Function >
void forRows (IndexType begin, IndexType end, Function &&function) const
 Method for parallel iteration over matrix rows from interval [begin, end) for constant instances.
 
template<typename Vector >
void getCompressedRowLengths (Vector &rowLengths) const
 Computes number of non-zeros in each row.
 
__cuda_callable__ DiagonalOffsetsView getDiagonalOffsets ()
 Returns a view for the vector of diagonal offsets.
 
__cuda_callable__ DiagonalOffsetsView::ConstViewType getDiagonalOffsets () const
 Returns a view for the vector of diagonal offsets.
 
__cuda_callable__ IndexType getDiagonalsCount () const
 Returns number of diagonals.
 
__cuda_callable__ RealType getElement (IndexType row, IndexType column) const
 Returns value of matrix element at position given by its row and column index.
 
__cuda_callable__ IndexerTypegetIndexer ()
 This method returns matrix elements indexer used by this matrix.
 
__cuda_callable__ const IndexerTypegetIndexer () const
 This method returns matrix elements indexer used by this matrix.
 
IndexType getNonzeroElementsCount () const override
 Returns number of non-zero matrix elements.
 
__cuda_callable__ RowView getRow (IndexType rowIdx)
 Non-constant getter of simple structure for accessing given matrix row.
 
__cuda_callable__ ConstRowView getRow (IndexType rowIdx) const
 Constant getter of simple structure for accessing given matrix row.
 
template<typename Vector >
void getRowCapacities (Vector &rowCapacities) const
 Compute capacities of all rows.
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_>
bool operator!= (const MultidiagonalMatrixBase< Real_, Device_, Index_, Organization_ > &matrix) const
 Comparison operator with another multidiagonal matrix.
 
MultidiagonalMatrixBaseoperator= (const MultidiagonalMatrixBase &)=delete
 Copy-assignment operator.
 
MultidiagonalMatrixBaseoperator= (MultidiagonalMatrixBase &&)=delete
 Move-assignment operator.
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_>
bool operator== (const MultidiagonalMatrixBase< Real_, Device_, Index_, Organization_ > &matrix) const
 Comparison operator with another multidiagonal matrix.
 
void print (std::ostream &str) const
 Method for printing the matrix to output stream.
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void reduceAllRows (Fetch &fetch, Reduce &reduce, Keep &keep, const FetchReal &identity) const
 Method for performing general reduction on all matrix rows for constant instances.
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void reduceRows (IndexType begin, IndexType end, Fetch &fetch, Reduce &reduce, Keep &keep, const FetchReal &identity) const
 Method for performing general reduction on matrix rows for constant instances.
 
template<typename Function >
void sequentialForAllRows (Function &function)
 This method calls sequentialForRows for all matrix rows.
 
template<typename Function >
void sequentialForAllRows (Function &function) const
 This method calls sequentialForRows for all matrix rows (for constant instances).
 
template<typename Function >
void sequentialForRows (IndexType begin, IndexType end, Function &function)
 Method for sequential iteration over all matrix rows for non-constant instances.
 
template<typename Function >
void sequentialForRows (IndexType begin, IndexType end, Function &function) const
 Method for sequential iteration over all matrix rows for constant instances.
 
__cuda_callable__ void setElement (IndexType row, IndexType column, const RealType &value)
 Sets element at given row and column to given value.
 
void setValue (const RealType &value)
 Set all matrix elements to given value.
 
template<typename InVector , typename OutVector >
void vectorProduct (const InVector &inVector, OutVector &outVector, RealType matrixMultiplicator=1.0, RealType outVectorMultiplicator=0.0, IndexType begin=0, IndexType end=0) const
 Computes product of matrix and vector.
 
- Public Member Functions inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
__cuda_callable__ MatrixBase ()=default
 Basic constructor with no parameters.
 
__cuda_callable__ MatrixBase (const MatrixBase &view)=default
 Shallow copy constructor.
 
__cuda_callable__ MatrixBase (IndexType rows, IndexType columns, ValuesViewType values)
 Constructor with matrix dimensions and matrix elements values.
 
__cuda_callable__ MatrixBase (MatrixBase &&view) noexcept=default
 Move constructor.
 
IndexType getAllocatedElementsCount () const
 Tells the number of allocated matrix elements.
 
__cuda_callable__ IndexType getColumns () const
 Returns number of matrix columns.
 
__cuda_callable__ IndexType getRows () const
 Returns number of matrix rows.
 
__cuda_callable__ ValuesViewTypegetValues ()
 Returns a reference to a vector with the matrix elements values.
 
__cuda_callable__ const ValuesViewTypegetValues () const
 Returns a constant reference to a vector with the matrix elements values.
 
bool operator!= (const Matrix &matrix) const
 Comparison operator with another arbitrary matrix view type.
 
bool operator!= (const MatrixT &matrix) const
 
__cuda_callable__ MatrixBaseoperator= (const MatrixBase &)=delete
 Copy-assignment operator.
 
__cuda_callable__ MatrixBaseoperator= (MatrixBase &&)=delete
 Move-assignment operator.
 
bool operator== (const Matrix &matrix) const
 Comparison operator with another arbitrary matrix view type.
 
bool operator== (const MatrixT &matrix) const
 

Additional Inherited Members

- Static Public Member Functions inherited from TNL::Matrices::MultidiagonalMatrixBase< Real, Device, Index, Organization >
static std::string getSerializationType ()
 Returns string with serialization type.
 
- Static Public Member Functions inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
static constexpr ElementsOrganization getOrganization ()
 Matrix elements organization getter.
 
static constexpr bool isBinary ()
 Test of binary matrix type.
 
static constexpr bool isMatrix ()
 Test of matrix type.
 
static constexpr bool isSymmetric ()
 Test of symmetric matrix type.
 
- Protected Member Functions inherited from TNL::Matrices::MultidiagonalMatrixBase< Real, Device, Index, Organization >
__cuda_callable__ void bind (typename Base::ValuesViewType values, DiagonalOffsetsView diagonalOffsets, HostDiagonalOffsetsView hostDiagonalOffsets, IndexerType indexer)
 Re-initializes the internal attributes of the base class.
 
- Protected Member Functions inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
__cuda_callable__ void bind (IndexType rows, IndexType columns, ValuesViewType values)
 Re-initializes the internal attributes of the base class.
 
- Protected Attributes inherited from TNL::Matrices::MultidiagonalMatrixBase< Real, Device, Index, Organization >
DiagonalOffsetsView diagonalOffsets
 
HostDiagonalOffsetsView hostDiagonalOffsets
 
IndexerType indexer
 
- Protected Attributes inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
IndexType columns
 
IndexType rows
 
ValuesViewType values
 

Detailed Description

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization()>
class TNL::Matrices::MultidiagonalMatrixView< Real, Device, Index, Organization >

Implementation of sparse multidiagonal matrix.

It serves as an accessor to SparseMatrix for example when passing the matrix to lambda functions. SparseMatrix view can be also created in CUDA kernels.

See MultidiagonalMatrix for more details.

Template Parameters
Realis a type of matrix elements.
Deviceis a device where the matrix is allocated.
Indexis a type for indexing of the matrix elements.
Organizationtells the ordering of matrix elements. It is either RowMajorOrder or ColumnMajorOrder.

Constructor & Destructor Documentation

◆ MultidiagonalMatrixView()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
__cuda_callable__ TNL::Matrices::MultidiagonalMatrixView< Real, Device, Index, Organization >::MultidiagonalMatrixView ( typename Base::ValuesViewType values,
typename Base::DiagonalOffsetsView diagonalOffsets,
typename Base::HostDiagonalOffsetsView hostDiagonalOffsets,
typename Base::IndexerType indexer )

Constructor with all necessary data and views.

Parameters
valuesis a vector view with matrix elements values
diagonalOffsetsis a vector view with diagonals offsets
hostDiagonalOffsetsis a vector view with a copy of diagonals offsets on the host
indexeris an indexer of matrix elements

Member Function Documentation

◆ bind()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
__cuda_callable__ void TNL::Matrices::MultidiagonalMatrixView< Real, Device, Index, Organization >::bind ( MultidiagonalMatrixView< Real, Device, Index, Organization > view)

Method for rebinding (reinitialization) using another multidiagonal matrix view.

Parameters
viewThe multidiagonal matrix view to be bound.

◆ getConstView()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
auto TNL::Matrices::MultidiagonalMatrixView< Real, Device, Index, Organization >::getConstView ( ) const

Returns a non-modifiable view of the multidiagonal matrix.

Returns
multidiagonal matrix view.

◆ getView()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
auto TNL::Matrices::MultidiagonalMatrixView< Real, Device, Index, Organization >::getView ( )

Returns a modifiable view of the multidiagonal matrix.

Returns
multidiagonal matrix view.

◆ operator=()

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization()>
MultidiagonalMatrixView & TNL::Matrices::MultidiagonalMatrixView< Real, Device, Index, Organization >::operator= ( const MultidiagonalMatrixView< Real, Device, Index, Organization > & )
delete

Copy-assignment operator.

It is a deleted function, because matrix assignment in general requires reallocation.

◆ save()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
void TNL::Matrices::MultidiagonalMatrixView< Real, Device, Index, Organization >::save ( File & file) const

Method for saving the matrix to a file.

Parameters
fileis the output file.

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