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

Implementation of sparse multidiagonal matrix. More...

#include <TNL/Matrices/MultidiagonalMatrix.h>

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

Public Types

using ConstViewType = MultidiagonalMatrixView< std::add_const_t< Real >, Device, Index, Organization >
 Matrix view type for constant instances.
 
using DiagonalOffsetsType = Containers::Vector< Index, Device, Index, IndexAllocator >
 
using HostDiagonalOffsetsType = Containers::Vector< Index, Devices::Host, Index >
 
using IndexAllocatorType = IndexAllocator
 The allocator for matrix elements offsets from the diagonal.
 
using RealAllocatorType = RealAllocator
 The allocator for matrix elements values.
 
template<typename _Real = Real, typename _Device = Device, typename _Index = Index, ElementsOrganization _Organization = Organization, typename _RealAllocator = typename Allocators::Default< _Device >::template Allocator< _Real >, typename _IndexAllocator = typename Allocators::Default< _Device >::template Allocator< _Index >>
using Self = MultidiagonalMatrix< _Real, _Device, _Index, _Organization, _RealAllocator, _IndexAllocator >
 Helper type for getting self type or its modifications.
 
using ValuesVectorType = Containers::Vector< Real, Device, Index, RealAllocator >
 Type of vector holding values of matrix elements.
 
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

 MultidiagonalMatrix ()=default
 Constructor with no parameters.
 
 MultidiagonalMatrix (const MultidiagonalMatrix &matrix)
 Copy constructor.
 
template<typename ListIndex , typename ListReal >
 MultidiagonalMatrix (Index columns, std::initializer_list< ListIndex > diagonalOffsets, const std::initializer_list< std::initializer_list< ListReal > > &data)
 Constructor with matrix dimensions, diagonals offsets and matrix elements.
 
 MultidiagonalMatrix (Index rows, Index columns)
 Constructor with matrix dimensions.
 
template<typename Vector >
 MultidiagonalMatrix (Index rows, Index columns, const Vector &diagonalOffsets)
 Constructor with matrix dimensions and matrix elements offsets.
 
template<typename ListIndex >
 MultidiagonalMatrix (Index rows, Index columns, std::initializer_list< ListIndex > diagonalOffsets)
 Constructor with matrix dimensions and diagonals offsets.
 
 MultidiagonalMatrix (MultidiagonalMatrix &&) noexcept=default
 Move constructor.
 
ConstViewType getConstView () const
 Returns a non-modifiable view of the multidiagonal matrix.
 
template<typename Real2 , typename Index2 >
void getTransposition (const MultidiagonalMatrix< Real2, Device, Index2 > &matrix, const Real &matrixMultiplicator=1.0)
 
ViewType getView ()
 Returns a modifiable view of the multidiagonal matrix.
 
void load (const String &fileName)
 Method for loading the matrix from the file with given filename.
 
void load (File &file) override
 Method for loading the matrix from a file.
 
MultidiagonalMatrixoperator= (const MultidiagonalMatrix &matrix)
 Assignment of exactly the same matrix type.
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_, typename RealAllocator_ , typename IndexAllocator_ >
MultidiagonalMatrixoperator= (const MultidiagonalMatrix< Real_, Device_, Index_, Organization_, RealAllocator_, IndexAllocator_ > &matrix)
 Assignment of another multidiagonal matrix.
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_, typename RealAllocator_ , typename IndexAllocator_ >
MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator > & operator= (const MultidiagonalMatrix< Real_, Device_, Index_, Organization_, RealAllocator_, IndexAllocator_ > &matrix)
 
void reset ()
 Resets the matrix to zero dimensions.
 
void save (const String &fileName) const
 Method for saving the matrix to the file with given filename.
 
void save (File &file) const override
 Method for saving the matrix to a file.
 
template<typename Vector >
void setDiagonalOffsets (const Vector &diagonalOffsets)
 Set the diagonals offsets by means of vector-like container.
 
template<typename ListIndex >
void setDiagonalOffsets (std::initializer_list< ListIndex > diagonalOffsets)
 Set the diagonals offsets by means of initializer list.
 
template<typename Vector >
void setDimensions (Index rows, Index columns, const Vector &diagonalOffsets)
 Set matrix dimensions and diagonals offsets.
 
template<typename ListReal >
void setElements (const std::initializer_list< std::initializer_list< ListReal > > &data)
 Set matrix elements from an initializer list.
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_, typename RealAllocator_ , typename IndexAllocator_ >
void setLike (const MultidiagonalMatrix< Real_, Device_, Index_, Organization_, RealAllocator_, IndexAllocator_ > &matrix)
 Setup the matrix dimensions and diagonals offsets based on another multidiagonal matrix.
 
template<typename RowCapacitiesVector >
void setRowCapacities (const RowCapacitiesVector &rowCapacities)
 This method is for compatibility with SparseMatrix.
 
- Public Member Functions inherited from TNL::Object
virtual ~Object ()=default
 Destructor.
 
virtual std::string getSerializationTypeVirtual () const
 
void load (const String &fileName)
 Method for restoring the object from a file.
 
void save (const String &fileName) const
 Method for saving the object to a file as a binary data.
 
- 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
 

Static Public Member Functions

static std::string getSerializationType ()
 Returns string with serialization type.
 
- Static Public Member Functions inherited from TNL::Object
static std::string getSerializationType ()
 Static serialization type getter.
 
- 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 Attributes

DiagonalOffsetsType diagonalOffsets
 
HostDiagonalOffsetsType hostDiagonalOffsets
 
ValuesVectorType values
 Vector containing the allocated matrix elements.
 
- 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
 

Additional Inherited Members

- 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.
 

Detailed Description

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
class TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >

Implementation of sparse multidiagonal matrix.

Use this matrix type for storing of matrices where the offsets of non-zero elements from the diagonal are the same in each row. Typically such matrices arise from discretization of partial differential equations on regular numerical grids. This is one example (dots represent zero matrix elements):

\[ \left( \begin{array}{ccccccc} 4 & -1 & . & -1 & . & . \\ -1 & 4 & -1 & . & -1 & . \\ . & -1 & 4 & -1 & . & -1 \\ -1 & . & -1 & 4 & -1 & . \\ . & -1 & . & -1 & 4 & -1 \\ . & . & -1 & . & -1 & 4 \end{array} \right) \]

In this matrix, the column indexes in each row \(i\) can be expressed as \(\{i-3, i-1, i, i+1, i+3\}\) (where the resulting index is non-negative and smaller than the number of matrix columns). Therefore the diagonals offsets are \(\{-3,-1,0,1,3\}\). Advantage is that we do not store the column indexes explicitly as it is in SparseMatrix. This can reduce significantly the memory requirements which also means better performance. See the following table for the storage requirements comparison between TNL::Matrices::MultidiagonalMatrix and TNL::Matrices::SparseMatrix.

Real Index SparseMatrix MultidiagonalMatrix Ratio
float 32-bit int 8 bytes per element 4 bytes per element 50%
double 32-bit int 12 bytes per element 8 bytes per element 75%
float 64-bit int 12 bytes per element 4 bytes per element 30%
double 64-bit int 16 bytes per element 8 bytes per element 50%
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.
RealAllocatoris allocator for the matrix elements.
IndexAllocatoris allocator for the matrix elements offsets.

Member Typedef Documentation

◆ ConstViewType

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
using TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::ConstViewType = MultidiagonalMatrixView< std::add_const_t< Real >, Device, Index, Organization >

Matrix view type for constant instances.

See MultidiagonalMatrixView.

◆ ViewType

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
using TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::ViewType = MultidiagonalMatrixView< Real, Device, Index, Organization >

Type of related matrix view.

See MultidiagonalMatrixView.

Constructor & Destructor Documentation

◆ MultidiagonalMatrix() [1/4]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::MultidiagonalMatrix ( Index rows,
Index columns )

Constructor with matrix dimensions.

Parameters
rowsis number of matrix rows.
columnsis number of matrix columns.

◆ MultidiagonalMatrix() [2/4]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Vector >
TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::MultidiagonalMatrix ( Index rows,
Index columns,
const Vector & diagonalOffsets )

Constructor with matrix dimensions and matrix elements offsets.

Parameters
rowsis number of matrix rows.
columnsis number of matrix columns.
diagonalOffsetsare offsets of subdiagonals from the main diagonal.
Example
#include <iostream>
#include <TNL/Algorithms/parallelFor.h>
#include <TNL/Containers/StaticArray.h>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
laplaceOperatorMatrix()
{
/***
* Set matrix representing approximation of the Laplace operator on regular
* grid using the finite difference method.
*/
const int gridSize( 4 );
const int matrixSize = gridSize * gridSize;
TNL::Containers::Vector< int, Device > offsets{ -gridSize, -1, 0, 1, gridSize };
TNL::Matrices::MultidiagonalMatrix< double, Device > matrix( matrixSize, matrixSize, offsets );
auto matrixView = matrix.getView();
auto f = [ = ] __cuda_callable__( const TNL::Containers::StaticArray< 2, int >& i ) mutable
{
const int elementIdx = i[ 1 ] * gridSize + i[ 0 ];
auto row = matrixView.getRow( elementIdx );
if( i[ 0 ] == 0 || i[ 1 ] == 0 || i[ 0 ] == gridSize - 1 || i[ 1 ] == gridSize - 1 )
row.setElement( 2, 1.0 ); // set matrix elements corresponding to boundary grid nodes
// and Dirichlet boundary conditions, i.e. 1 on the main diagonal
// which is the third one
else {
row.setElement( 0, -1.0 ); // set matrix elements corresponding to inner grid nodes, i.e.
row.setElement( 1, -1.0 ); // 4 on the main diagonal (the third one) and -1 to the other
row.setElement( 2, 4.0 ); // sub-diagonals
row.setElement( 3, -1.0 );
row.setElement( 4, -1.0 );
}
};
TNL::Containers::StaticArray< 2, int > end = { gridSize, gridSize };
TNL::Algorithms::parallelFor< Device >( begin, end, f );
std::cout << "Laplace operator matrix: " << std::endl << matrix << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Creating Laplace operator matrix on CPU ... " << std::endl;
laplaceOperatorMatrix< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating Laplace operator matrix on CUDA GPU ... " << std::endl;
laplaceOperatorMatrix< TNL::Devices::Cuda >();
#endif
}
#define __cuda_callable__
Definition Macros.h:49
Array with constant size.
Definition StaticArray.h:20
Vector extends Array with algebraic operations.
Definition Vector.h:36
Implementation of sparse multidiagonal matrix.
Definition MultidiagonalMatrix.h:64
T end(T... args)
T endl(T... args)
Output
Creating Laplace operator matrix on CPU ...
Laplace operator matrix:
Row: 0 -> 0:1
Row: 1 -> 1:1
Row: 2 -> 2:1
Row: 3 -> 3:1
Row: 4 -> 4:1
Row: 5 -> 1:-1 4:-1 5:4 6:-1 9:-1
Row: 6 -> 2:-1 5:-1 6:4 7:-1 10:-1
Row: 7 -> 7:1
Row: 8 -> 8:1
Row: 9 -> 5:-1 8:-1 9:4 10:-1 13:-1
Row: 10 -> 6:-1 9:-1 10:4 11:-1 14:-1
Row: 11 -> 11:1
Row: 12 -> 12:1
Row: 13 -> 13:1
Row: 14 -> 14:1
Row: 15 -> 15:1
Creating Laplace operator matrix on CUDA GPU ...
Laplace operator matrix:
Row: 0 -> 0:1
Row: 1 -> 1:1
Row: 2 -> 2:1
Row: 3 -> 3:1
Row: 4 -> 4:1
Row: 5 -> 1:-1 4:-1 5:4 6:-1 9:-1
Row: 6 -> 2:-1 5:-1 6:4 7:-1 10:-1
Row: 7 -> 7:1
Row: 8 -> 8:1
Row: 9 -> 5:-1 8:-1 9:4 10:-1 13:-1
Row: 10 -> 6:-1 9:-1 10:4 11:-1 14:-1
Row: 11 -> 11:1
Row: 12 -> 12:1
Row: 13 -> 13:1
Row: 14 -> 14:1
Row: 15 -> 15:1

◆ MultidiagonalMatrix() [3/4]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename ListIndex >
TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::MultidiagonalMatrix ( Index rows,
Index columns,
std::initializer_list< ListIndex > diagonalOffsets )

Constructor with matrix dimensions and diagonals offsets.

Parameters
rowsis number of matrix rows.
columnsis number of matrix columns.
diagonalOffsetsare offsets of sub-diagonals from the main diagonal.
Example
#include <iostream>
#include <TNL/Algorithms/parallelFor.h>
#include <TNL/Containers/StaticArray.h>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
laplaceOperatorMatrix()
{
/***
* Set matrix representing approximation of the Laplace operator on regular
* grid using the finite difference method.
*/
const int gridSize( 4 );
const int matrixSize = gridSize * gridSize;
TNL::Matrices::MultidiagonalMatrix< double, Device > matrix( matrixSize, matrixSize, { -gridSize, -1, 0, 1, gridSize } );
auto matrixView = matrix.getView();
auto f = [ = ] __cuda_callable__( const TNL::Containers::StaticArray< 2, int >& i ) mutable
{
const int elementIdx = i[ 0 ] * gridSize + i[ 1 ];
auto row = matrixView.getRow( elementIdx );
if( i[ 0 ] == 0 || i[ 1 ] == 0 || i[ 0 ] == gridSize - 1 || i[ 1 ] == gridSize - 1 )
row.setElement( 2, 1.0 ); // set matrix elements corresponding to boundary grid nodes
// and Dirichlet boundary conditions, i.e. 1 on the main diagonal
// which is the third one
else {
row.setElement( 0, -1.0 ); // set matrix elements corresponding to inner grid nodes, i.e.
row.setElement( 1, -1.0 ); // 4 on the main diagonal (the third one) and -1 to the other
row.setElement( 2, 4.0 ); // sub-diagonals
row.setElement( 3, -1.0 );
row.setElement( 4, -1.0 );
}
};
TNL::Containers::StaticArray< 2, int > end = { gridSize, gridSize };
TNL::Algorithms::parallelFor< Device >( begin, end, f );
std::cout << "Laplace operator matrix: " << std::endl << matrix << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Creating Laplace operator matrix on CPU ... " << std::endl;
laplaceOperatorMatrix< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating Laplace operator matrix on CUDA GPU ... " << std::endl;
laplaceOperatorMatrix< TNL::Devices::Cuda >();
#endif
}
__cuda_callable__ RowView getRow(IndexType rowIdx)
Non-constant getter of simple structure for accessing given matrix row.
Definition MultidiagonalMatrixBase.hpp:145
__cuda_callable__ void setElement(IndexType localIdx, const RealType &value)
Changes value of matrix element on given subdiagonal.
Definition MultidiagonalMatrixRowView.hpp:67
ViewType getView()
Returns a modifiable view of the multidiagonal matrix.
Definition MultidiagonalMatrix.hpp:71
Output
Creating Laplace operator matrix on CPU ...
Laplace operator matrix:
Row: 0 -> 0:1
Row: 1 -> 1:1
Row: 2 -> 2:1
Row: 3 -> 3:1
Row: 4 -> 4:1
Row: 5 -> 1:-1 4:-1 5:4 6:-1 9:-1
Row: 6 -> 2:-1 5:-1 6:4 7:-1 10:-1
Row: 7 -> 7:1
Row: 8 -> 8:1
Row: 9 -> 5:-1 8:-1 9:4 10:-1 13:-1
Row: 10 -> 6:-1 9:-1 10:4 11:-1 14:-1
Row: 11 -> 11:1
Row: 12 -> 12:1
Row: 13 -> 13:1
Row: 14 -> 14:1
Row: 15 -> 15:1
Creating Laplace operator matrix on CUDA GPU ...
Laplace operator matrix:
Row: 0 -> 0:1
Row: 1 -> 1:1
Row: 2 -> 2:1
Row: 3 -> 3:1
Row: 4 -> 4:1
Row: 5 -> 1:-1 4:-1 5:4 6:-1 9:-1
Row: 6 -> 2:-1 5:-1 6:4 7:-1 10:-1
Row: 7 -> 7:1
Row: 8 -> 8:1
Row: 9 -> 5:-1 8:-1 9:4 10:-1 13:-1
Row: 10 -> 6:-1 9:-1 10:4 11:-1 14:-1
Row: 11 -> 11:1
Row: 12 -> 12:1
Row: 13 -> 13:1
Row: 14 -> 14:1
Row: 15 -> 15:1

◆ MultidiagonalMatrix() [4/4]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename ListIndex , typename ListReal >
TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::MultidiagonalMatrix ( Index columns,
std::initializer_list< ListIndex > diagonalOffsets,
const std::initializer_list< std::initializer_list< ListReal > > & data )

Constructor with matrix dimensions, diagonals offsets and matrix elements.

The number of matrix rows is deduced from the size of the initializer list data.

Template Parameters
ListIndexis type used in the initializer list defining matrix diagonals offsets.
ListRealis type used in the initializer list defining matrix elements values.
Parameters
columnsis number of matrix columns.
diagonalOffsetsare offsets of sub-diagonals from the main diagonal.
datais initializer list holding matrix elements. The size of the outer list defines the number of matrix rows. Each inner list defines values of each sub-diagonal and so its size should be lower or equal to the size of diagonalOffsets. Values of sub-diagonals which do not fit to given row are omitted.
Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
createMultidiagonalMatrix()
{
const int matrixSize = 6;
/*
* Setup the following matrix (dots represent zeros):
*
* / 4 -1 . -1 . . \
* | -1 4 -1 . -1 . |
* | . -1 4 -1 . -1 |
* | -1 . -1 4 -1 . |
* | . -1 . -1 4 -1 |
* \ . . 1 . -1 4 /
*
* The diagonals offsets are { -3, -1, 0, 1, 3 }.
*/
matrixSize,
{ -3, -1, 0, 1, 3 },
{
/*
* To set the matrix elements we first extend the diagonals to their full
* lengths even outside the matrix (dots represent zeros and zeros are
* artificial zeros used for memory alignment):
*
* 0 . 0 / 4 -1 . -1 . . \ -> { 0, 0, 4, -1, -1 }
* . 0 . | -1 4 -1 . -1 . | . -> { 0, -1, 4, -1, -1 }
* . . 0 | . -1 4 -1 . -1 | . . -> { 0, -1, 4, -1, -1 }
* . . | -1 . -1 4 -1 . | 0 . . -> { -1, -1, 4, -1, 0 }
* . | . -1 . -1 4 -1 | . 0 . . -> { -1, -1, 4, -1, 0 }
* \ . . 1 . -1 4 / 0 . 0 . . -> { -1, -1, 4, 0, 0 }
*/
// clang-format off
{ 0, 0, 4, -1, -1 },
{ 0, -1, 4, -1, -1 },
{ 0, -1, 4, -1, -1 },
{ -1, -1, 4, -1, 0 },
{ -1, -1, 4, -1, 0 },
{ -1, -1, 4, 0, 0 }
// clang-format on
} );
std::cout << "The matrix reads as: " << std::endl << matrix << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Create multidiagonal matrix on CPU ... " << std::endl;
createMultidiagonalMatrix< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating multidiagonal matrix on CUDA GPU ... " << std::endl;
createMultidiagonalMatrix< TNL::Devices::Cuda >();
#endif
}
Output
Create multidiagonal matrix on CPU ...
The matrix reads as:
Row: 0 -> 0:4 1:-1 3:-1
Row: 1 -> 0:-1 1:4 2:-1 4:-1
Row: 2 -> 1:-1 2:4 3:-1 5:-1
Row: 3 -> 0:-1 2:-1 3:4 4:-1
Row: 4 -> 1:-1 3:-1 4:4 5:-1
Row: 5 -> 2:-1 4:-1 5:4
Creating multidiagonal matrix on CUDA GPU ...
The matrix reads as:
Row: 0 -> 0:4 1:-1 3:-1
Row: 1 -> 0:-1 1:4 2:-1 4:-1
Row: 2 -> 1:-1 2:4 3:-1 5:-1
Row: 3 -> 0:-1 2:-1 3:4 4:-1
Row: 4 -> 1:-1 3:-1 4:4 5:-1
Row: 5 -> 2:-1 4:-1 5:4

Member Function Documentation

◆ getConstView()

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

Returns a non-modifiable view of the multidiagonal matrix.

See MultidiagonalMatrixView.

Returns
multidiagonal matrix view.

◆ getSerializationType()

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
std::string TNL::Matrices::MultidiagonalMatrixBase< Real, Device, Index, Organization >::getSerializationType ( )
static

Returns string with serialization type.

The string has a form Matrices::MultidiagonalMatrix< RealType, [any_device], IndexType, Organization, [any_allocator], [any_allocator] >.

See MultidiagonalMatrix::getSerializationType.

Returns
String with the serialization type.

◆ getView()

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

Returns a modifiable view of the multidiagonal matrix.

See MultidiagonalMatrixView.

Returns
multidiagonal matrix view.

◆ load() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::load ( const String & fileName)

Method for loading the matrix from the file with given filename.

Parameters
fileNameis name of the file.

◆ load() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::load ( File & file)
overridevirtual

Method for loading the matrix from a file.

Parameters
fileis the input file.

Reimplemented from TNL::Object.

◆ operator=() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator > & TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::operator= ( const MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator > & matrix)

Assignment of exactly the same matrix type.

Parameters
matrixis input matrix for the assignment.
Returns
reference to this matrix.

◆ operator=() [2/2]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_, typename RealAllocator_ , typename IndexAllocator_ >
MultidiagonalMatrix & TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::operator= ( const MultidiagonalMatrix< Real_, Device_, Index_, Organization_, RealAllocator_, IndexAllocator_ > & matrix)

Assignment of another multidiagonal matrix.

Parameters
matrixis input matrix for the assignment.
Returns
reference to this matrix.

◆ save() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::save ( const String & fileName) const

Method for saving the matrix to the file with given filename.

Parameters
fileNameis name of the file.

◆ save() [2/2]

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

Method for saving the matrix to a file.

Parameters
fileis the output file.

Reimplemented from TNL::Object.

◆ setDiagonalOffsets() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Vector >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::setDiagonalOffsets ( const Vector & diagonalOffsets)

Set the diagonals offsets by means of vector-like container.

This method deletes current matrix elements.

Template Parameters
Vectoris a type of vector-like container holding the subdiagonals offsets.
Parameters
diagonalOffsetsis a vector-like container holding the subdiagonals offsets.

◆ setDiagonalOffsets() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename ListIndex >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::setDiagonalOffsets ( std::initializer_list< ListIndex > diagonalOffsets)

Set the diagonals offsets by means of initializer list.

This method deletes current matrix elements.

Template Parameters
ListIndexis type of indexes used for the subdiagonals offsets definition.
Parameters
diagonalOffsetsis a initializer list with subdiagonals offsets.

◆ setDimensions()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Vector >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::setDimensions ( Index rows,
Index columns,
const Vector & diagonalOffsets )

Set matrix dimensions and diagonals offsets.

Template Parameters
Vectoris type of vector like container holding the diagonals offsets.
Parameters
rowsis number of matrix rows.
columnsis number of matrix columns.
diagonalOffsetsis vector with diagonals offsets.

◆ setElements()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename ListReal >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::setElements ( const std::initializer_list< std::initializer_list< ListReal > > & data)

Set matrix elements from an initializer list.

Template Parameters
ListRealis data type of the initializer list.
Parameters
datais initializer list holding matrix elements. The size of the outer list defines the number of matrix rows. Each inner list defines values of each sub-diagonal and so its size should be lower or equal to the size of diagonalsOffsets. Values of sub-diagonals which do not fit to given row are omitted.
Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
createMultidiagonalMatrix()
{
const int matrixSize = 6;
/*
* Setup the following matrix (dots represent zeros):
*
* / 4 -1 . -1 . . \
* | -1 4 -1 . -1 . |
* | . -1 4 -1 . -1 |
* | -1 . -1 4 -1 . |
* | . -1 . -1 4 -1 |
* \ . . 1 . -1 4 /
*
* The diagonals offsets are { -3, -1, 0, 1, 3 }.
*/
TNL::Matrices::MultidiagonalMatrix< double, Device > matrix( matrixSize, // number of matrix rows
matrixSize, // number of matrix columns
{ -3, -1, 0, 1, 3 } ); // matrix diagonals offsets
/*
* To set the matrix elements we first extend the diagonals to their full
* lengths even outside the matrix (dots represent zeros and zeros are
* artificial zeros used for memory alignment):
*
* 0 . 0 / 4 -1 . -1 . . \ -> { 0, 0, 4, -1, -1 }
* . 0 . | -1 4 -1 . -1 . | . -> { 0, -1, 4, -1, -1 }
* . . 0 | . -1 4 -1 . -1 | . . -> { 0, -1, 4, -1, -1 }
* . . | -1 . -1 4 -1 . | 0 . . -> { -1, -1, 4, -1, 0 }
* . | . -1 . -1 4 -1 | . 0 . . -> { -1, -1, 4, -1, 0 }
* \ . . 1 . -1 4 / 0 . 0 . . -> { -1, -1, 4, 0, 0 }
*/
matrix.setElements( {
// clang-format off
{ 0, 0, 4, -1, -1 },
{ 0, -1, 4, -1, -1 },
{ 0, -1, 4, -1, -1 },
{ -1, -1, 4, -1, 0 },
{ -1, -1, 4, -1, 0 },
{ -1, -1, 4, 0, 0 }
// clang-format on
} );
std::cout << "The matrix reads as: " << std::endl << matrix << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Create multidiagonal matrix on CPU ... " << std::endl;
createMultidiagonalMatrix< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating multidiagonal matrix on CUDA GPU ... " << std::endl;
createMultidiagonalMatrix< TNL::Devices::Cuda >();
#endif
}
void setElements(const std::initializer_list< std::initializer_list< ListReal > > &data)
Set matrix elements from an initializer list.
Definition MultidiagonalMatrix.hpp:203
Output
Create multidiagonal matrix on CPU ...
The matrix reads as:
Row: 0 -> 0:4 1:-1 3:-1
Row: 1 -> 0:-1 1:4 2:-1 4:-1
Row: 2 -> 1:-1 2:4 3:-1 5:-1
Row: 3 -> 0:-1 2:-1 3:4 4:-1
Row: 4 -> 1:-1 3:-1 4:4 5:-1
Row: 5 -> 2:-1 4:-1 5:4
Creating multidiagonal matrix on CUDA GPU ...
The matrix reads as:
Row: 0 -> 0:4 1:-1 3:-1
Row: 1 -> 0:-1 1:4 2:-1 4:-1
Row: 2 -> 1:-1 2:4 3:-1 5:-1
Row: 3 -> 0:-1 2:-1 3:4 4:-1
Row: 4 -> 1:-1 3:-1 4:4 5:-1
Row: 5 -> 2:-1 4:-1 5:4

◆ setLike()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_, typename RealAllocator_ , typename IndexAllocator_ >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::setLike ( const MultidiagonalMatrix< Real_, Device_, Index_, Organization_, RealAllocator_, IndexAllocator_ > & matrix)

Setup the matrix dimensions and diagonals offsets based on another multidiagonal matrix.

Template Parameters
Real_is Real type of the source matrix.
Device_is Device type of the source matrix.
Index_is Index type of the source matrix.
Organization_is Organization of the source matrix.
RealAllocator_is RealAllocator of the source matrix.
IndexAllocator_is IndexAllocator of the source matrix.
Parameters
matrixis the source matrix.

◆ setRowCapacities()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename RowCapacitiesVector >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::setRowCapacities ( const RowCapacitiesVector & rowCapacities)

This method is for compatibility with SparseMatrix.

It checks if the number of matrix diagonals is compatible with required number of non-zero matrix elements in each row. If not exception is thrown.

Template Parameters
RowCapacitiesVectoris vector-like container type for holding required row capacities.
Parameters
rowCapacitiesis vector-like container holding required row capacities.

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