Template Numerical Library version\ main:8b8c8226
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 BaseType = Matrix< Real, Device, Index, RealAllocator >
 
using ConstRowView = typename ViewType::ConstRowView
 Type for accessing constant matrix rows.
 
using ConstViewType = MultidiagonalMatrixView< std::add_const_t< Real >, Device, Index, Organization >
 Matrix view type for constant instances. More...
 
using DeviceType = Device
 The device where the matrix is allocated.
 
using DiagonalsOffsetsType = Containers::Vector< Index, Device, Index, IndexAllocator >
 
using DiagonalsOffsetsView = typename DiagonalsOffsetsType::ViewType
 
using HostDiagonalsOffsetsType = Containers::Vector< Index, Devices::Host, Index >
 
using HostDiagonalsOffsetsView = typename HostDiagonalsOffsetsType::ViewType
 
using IndexAllocatorType = IndexAllocator
 The allocator for matrix elements offsets from the diagonal.
 
using IndexerType = details::MultidiagonalMatrixIndexer< Index, Organization==Algorithms::Segments::RowMajorOrder >
 
using IndexType = Index
 The type used for matrix elements indexing.
 
using RealAllocatorType = RealAllocator
 The allocator for matrix elements values.
 
using RealType = Real
 The type of matrix elements.
 
using RowView = typename ViewType::RowView
 Type for accessing matrix rows.
 
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 ValuesType = typename BaseType::ValuesType
 
using ValuesView = typename ValuesType::ViewType
 
using ViewType = MultidiagonalMatrixView< Real, Device, Index, Organization >
 Type of related matrix view. More...
 
- Public Types inherited from TNL::Matrices::Matrix< Real, Device, Index, RealAllocator >
using ConstRowsCapacitiesView = typename RowsCapacitiesView::ConstViewType
 
using ConstValuesType = Containers::Vector< std::add_const_t< Real >, Device, Index, RealAllocator >
 Type of constant vector holding values of matrix elements.
 
using ConstValuesView = typename ViewType::ConstValuesView
 Type of constant vector view holding values of matrix elements.
 
using ConstViewType = typename MatrixView< Real, Device, Index >::ConstViewType
 Type of base matrix view for constant instances. More...
 
using DeviceType = Device
 The device where the matrix is allocated.
 
using IndexType = Index
 The type used for matrix elements indexing.
 
using RealAllocatorType = RealAllocator
 
using RealType = Real
 The type of matrix elements.
 
using RowsCapacitiesType = Containers::Vector< Index, Device, Index >
 
using RowsCapacitiesView = Containers::VectorView< Index, Device, Index >
 
using ValuesType = Containers::Vector< Real, Device, Index, RealAllocator >
 Type of vector holding values of matrix elements.
 
using ValuesView = typename ViewType::ValuesView
 Type of vector view holding values of matrix elements.
 
using ViewType = MatrixView< Real, Device, Index >
 Type of base matrix view. More...
 

Public Member Functions

 MultidiagonalMatrix ()
 Constructor with no parameters.
 
 MultidiagonalMatrix (const MultidiagonalMatrix &matrix)=default
 Copy constructor. More...
 
template<typename ListIndex , typename ListReal >
 MultidiagonalMatrix (IndexType columns, std::initializer_list< ListIndex > diagonalsOffsets, const std::initializer_list< std::initializer_list< ListReal > > &data)
 Constructor with matrix dimensions, diagonals offsets and matrix elements. More...
 
 MultidiagonalMatrix (IndexType rows, IndexType columns)
 Constructor with matrix dimensions. More...
 
template<typename Vector >
 MultidiagonalMatrix (IndexType rows, IndexType columns, const Vector &diagonalsOffsets)
 Constructor with matrix dimensions and matrix elements offsets. More...
 
template<typename ListIndex >
 MultidiagonalMatrix (IndexType rows, IndexType columns, std::initializer_list< ListIndex > diagonalsOffsets)
 Constructor with matrix dimensions and diagonals offsets. More...
 
 MultidiagonalMatrix (MultidiagonalMatrix &&matrix) noexcept=default
 Move constructor. More...
 
__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. More...
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_, typename RealAllocator_ >
void addMatrix (const MultidiagonalMatrix< Real_, Device_, Index_, Organization_, RealAllocator_ > &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. More...
 
template<typename Function >
void forAllElements (Function &function) const
 This method calls forElements for all matrix rows (for constant instances). More...
 
template<typename Function >
void forAllRows (Function &&function)
 Method for parallel iteration over all matrix rows. More...
 
template<typename Function >
void forAllRows (Function &&function) const
 Method for parallel iteration over all matrix rows for constant instances. More...
 
template<typename Function >
void forElements (IndexType begin, IndexType end, Function &function)
 Method for iteration over matrix rows for non-constant instances. More...
 
template<typename Function >
void forElements (IndexType begin, IndexType end, Function &function) const
 Method for iteration over matrix rows for constant instances. More...
 
template<typename Function >
void forRows (IndexType begin, IndexType end, Function &&function)
 Method for parallel iteration over matrix rows from interval [ begin, end). More...
 
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. More...
 
template<typename Vector >
void getCompressedRowLengths (Vector &rowLengths) const
 Computes number of non-zeros in each row. More...
 
ConstViewType getConstView () const
 Returns a non-modifiable view of the multidiagonal matrix. More...
 
IndexType getDiagonalsCount () const
 Returns number of diagonals. More...
 
const DiagonalsOffsetsTypegetDiagonalsOffsets () const
 Returns vector with diagonals offsets. More...
 
__cuda_callable__ RealType getElement (IndexType row, IndexType column) const
 Returns value of matrix element at position given by its row and column index. More...
 
IndexerTypegetIndexer ()
 This method returns matrix elements indexer used by this matrix. More...
 
const IndexerTypegetIndexer () const
 This method returns matrix elements indexer used by this matrix. More...
 
IndexType getNonzeroElementsCount () const override
 Returns number of non-zero matrix elements. More...
 
__cuda_callable__ IndexType getPaddingIndex () const
 Returns padding index denoting padding zero elements. More...
 
__cuda_callable__ RowView getRow (IndexType rowIdx)
 Non-constant getter of simple structure for accessing given matrix row. More...
 
__cuda_callable__ ConstRowView getRow (IndexType rowIdx) const
 Constant getter of simple structure for accessing given matrix row. More...
 
template<typename Vector >
void getRowCapacities (Vector &rowCapacities) const
 Compute capacities of all rows. More...
 
IndexType getRowLength (IndexType row) const
 
std::string getSerializationTypeVirtual () const override
 Returns string with serialization type. More...
 
template<typename Real2 , typename Index2 >
void getTransposition (const MultidiagonalMatrix< Real2, Device, Index2 > &matrix, const RealType &matrixMultiplicator=1.0)
 
ViewType getView ()
 Returns a modifiable view of the multidiagonal matrix. More...
 
void load (const String &fileName)
 Method for loading the matrix from the file with given filename. More...
 
void load (File &file) override
 Method for loading the matrix from a file. More...
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_, typename RealAllocator_ , typename IndexAllocator_ >
bool operator!= (const MultidiagonalMatrix< Real_, Device_, Index_, Organization_, RealAllocator_, IndexAllocator_ > &matrix) const
 Comparison operator with another multidiagonal matrix. More...
 
MultidiagonalMatrixoperator= (const MultidiagonalMatrix &matrix)
 Assignment of exactly the same matrix type. More...
 
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. More...
 
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)
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_, typename RealAllocator_ , typename IndexAllocator_ >
bool operator== (const MultidiagonalMatrix< Real_, Device_, Index_, Organization_, RealAllocator_, IndexAllocator_ > &matrix) const
 Comparison operator with another multidiagonal matrix. More...
 
void print (std::ostream &str) const override
 Method for printing the matrix to output stream. More...
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void reduceAllRows (Fetch &fetch, Reduce &reduce, Keep &keep, const FetchReal &identity)
 Method for performing general reduction on all matrix rows. More...
 
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. More...
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void reduceRows (IndexType begin, IndexType end, Fetch &fetch, Reduce &reduce, Keep &keep, const FetchReal &identity)
 Method for performing general reduction on matrix rows. More...
 
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. More...
 
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. More...
 
void save (File &file) const override
 Method for saving the matrix to a file. More...
 
template<typename Function >
void sequentialForAllRows (Function &function)
 This method calls sequentialForRows for all matrix rows. More...
 
template<typename Function >
void sequentialForAllRows (Function &function) const
 This method calls sequentialForRows for all matrix rows (for constant instances). More...
 
template<typename Function >
void sequentialForRows (IndexType begin, IndexType end, Function &function)
 Method for sequential iteration over all matrix rows for non-constant instances. More...
 
template<typename Function >
void sequentialForRows (IndexType begin, IndexType end, Function &function) const
 Method for sequential iteration over all matrix rows for constant instances. More...
 
template<typename Vector >
void setDiagonalsOffsets (const Vector &diagonalsOffsets)
 Set the diagonals offsets by means of vector-like container. More...
 
template<typename ListIndex >
void setDiagonalsOffsets (std::initializer_list< ListIndex > diagonalsOffsets)
 Set the diagonals offsets by means of initializer list. More...
 
template<typename Vector >
void setDimensions (IndexType rows, IndexType columns, const Vector &diagonalsOffsets)
 Set matrix dimensions and diagonals offsets. More...
 
__cuda_callable__ void setElement (IndexType row, IndexType column, const RealType &value)
 Sets element at given row and column to given value. More...
 
template<typename ListReal >
void setElements (const std::initializer_list< std::initializer_list< ListReal > > &data)
 Set matrix elements from an initializer list. More...
 
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. More...
 
template<typename RowCapacitiesVector >
void setRowCapacities (const RowCapacitiesVector &rowCapacities)
 This method is for compatibility with SparseMatrix. More...
 
void setValue (const RealType &value)
 Set all matrix elements to given value. More...
 
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. More...
 
- Public Member Functions inherited from TNL::Matrices::Matrix< Real, Device, Index, RealAllocator >
 Matrix (const RealAllocatorType &allocator=RealAllocatorType())
 Construct a new Matrix object possibly with user defined allocator of the matrix values. More...
 
 Matrix (IndexType rows, IndexType columns, const RealAllocatorType &allocator=RealAllocatorType())
 Construct a new Matrix object with given dimensions and possibly user defined allocator of the matrix values. More...
 
IndexType getAllocatedElementsCount () const
 Tells the number of allocated matrix elements. More...
 
__cuda_callable__ IndexType getColumns () const
 Returns number of matrix columns. More...
 
virtual IndexType getNonzeroElementsCount () const
 Computes a current number of nonzero matrix elements. More...
 
__cuda_callable__ IndexType getRows () const
 Returns number of matrix rows. More...
 
ValuesTypegetValues ()
 Returns a reference to a vector with the matrix elements values. More...
 
const ValuesTypegetValues () const
 Returns a constant reference to a vector with the matrix elements values. More...
 
void load (File &file) override
 Method for loading the matrix from a file. More...
 
template<typename Matrix >
bool operator!= (const Matrix &matrix) const
 Comparison operator with another arbitrary matrix type. More...
 
template<typename MatrixT >
bool operator!= (const MatrixT &matrix) const
 
template<typename Matrix >
bool operator== (const Matrix &matrix) const
 Comparison operator with another arbitrary matrix type. More...
 
template<typename MatrixT >
bool operator== (const MatrixT &matrix) const
 
virtual void print (std::ostream &str) const
 Method for printing the matrix to output stream. More...
 
void reset ()
 Reset the matrix. More...
 
void save (File &file) const override
 Method for saving the matrix to a file. More...
 
virtual void setDimensions (IndexType rows, IndexType columns)
 Method for setting or changing of the matrix dimensions. More...
 
template<typename Matrix_ >
void setLike (const Matrix_ &matrix)
 Set the matrix dimensions to be equal to those of the input matrix. More...
 
- Public Member Functions inherited from TNL::Object
virtual ~Object ()=default
 Destructor. More...
 
virtual std::string getSerializationTypeVirtual () const
 
void load (const String &fileName)
 Method for restoring the object from a file. More...
 
virtual void load (File &file)
 Method for restoring the object from a file. More...
 
void save (const String &fileName) const
 Method for saving the object to a file as a binary data. More...
 
virtual void save (File &file) const
 Method for saving the object to a file as a binary data. More...
 

Static Public Member Functions

static constexpr ElementsOrganization getOrganization ()
 Elements organization getter.
 
static std::string getSerializationType ()
 Returns string with serialization type. More...
 
static constexpr bool isSymmetric ()
 This is only for compatibility with sparse matrices. More...
 
- Static Public Member Functions inherited from TNL::Object
static std::string getSerializationType ()
 Static serialization type getter. More...
 

Protected Attributes

DiagonalsOffsetsType diagonalsOffsets
 
HostDiagonalsOffsetsType hostDiagonalsOffsets
 
IndexerType indexer
 
ViewType view
 
- Protected Attributes inherited from TNL::Matrices::Matrix< Real, Device, Index, RealAllocator >
IndexType columns
 
IndexType rows
 
ValuesType values
 Array containing the allocated matrix elements.
 

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/6]

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 ( IndexType  rows,
IndexType  columns 
)

Constructor with matrix dimensions.

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

◆ MultidiagonalMatrix() [2/6]

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 ( IndexType  rows,
IndexType  columns,
const Vector &  diagonalsOffsets 
)

Constructor with matrix dimensions and matrix elements offsets.

Parameters
rowsis number of matrix rows.
columnsis number of matrix columns.
diagonalsOffsetsare 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: CudaCallable.h:22
Array with constant size.
Definition: StaticArray.h:24
Vector extends Array with algebraic operations.
Definition: Vector.h:40
Implementation of sparse multidiagonal matrix.
Definition: MultidiagonalMatrix.h:71
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/6]

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 ( IndexType  rows,
IndexType  columns,
std::initializer_list< ListIndex >  diagonalsOffsets 
)

Constructor with matrix dimensions and diagonals offsets.

Parameters
rowsis number of matrix rows.
columnsis number of matrix columns.
diagonalsOffsetsare 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: MultidiagonalMatrix.hpp:422
ViewType getView()
Returns a modifiable view of the multidiagonal matrix.
Definition: MultidiagonalMatrix.hpp:83
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/6]

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 ( IndexType  columns,
std::initializer_list< ListIndex >  diagonalsOffsets,
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.
diagonalsOffsetsare 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 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 }.
*/
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 }
*
*/
{ 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 }
} );
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

◆ MultidiagonalMatrix() [5/6]

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 ( const MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator > &  matrix)
default

Copy constructor.

Parameters
matrixis an input matrix.

◆ MultidiagonalMatrix() [6/6]

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 ( MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator > &&  matrix)
defaultnoexcept

Move constructor.

Parameters
matrixis an input matrix.

Member Function Documentation

◆ addElement()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
__cuda_callable__ void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::addElement ( IndexType  row,
IndexType  column,
const RealType value,
const RealType thisElementMultiplicator = 1.0 
)

Add element at given row and column to given value.

This method can be called from the host system (CPU) no matter where the matrix is allocated. If the matrix is allocated on GPU this method can be called even from device kernels. If the matrix is allocated in GPU device this method is called from CPU, it transfers values of each matrix element separately and so the performance is very low. For higher performance see. MultidiagonalMatrix::getRow or MultidiagonalMatrix::forElements and MultidiagonalMatrix::forAllElements. The call may fail if the matrix row capacity is exhausted.

Parameters
rowis row index of the element.
columnis columns index of the element.
valueis the value the element will be set to.
thisElementMultiplicatoris multiplicator the original matrix element value is multiplied by before addition of given value.
Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void addElements()
{
const int matrixSize( 5 );
matrixSize, // number of rows
matrixSize, // number of columns
{ -1, 0, 1 } ); // diagonals offsets
for( int i = 0; i < matrixSize; i++ )
matrix.setElement( i, i, i );
std::cout << "Initial matrix is: " << std::endl << matrix << std::endl;
for( int i = 0; i < matrixSize; i++ )
{
if( i > 0 )
matrix.addElement( i, i - 1, 1.0, 5.0 );
matrix.addElement( i, i, 1.0, 5.0 );
if( i < matrixSize - 1 )
matrix.addElement( i, i + 1, 1.0, 5.0 );
}
std::cout << "Matrix after addition is: " << std::endl << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Add elements on host:" << std::endl;
addElements< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Add elements on CUDA device:" << std::endl;
addElements< TNL::Devices::Cuda >();
#endif
}
__cuda_callable__ void setElement(IndexType row, IndexType column, const RealType &value)
Sets element at given row and column to given value.
Definition: MultidiagonalMatrix.hpp:435
Output
Add elements on host:
Initial matrix is:
Row: 0 ->
Row: 1 -> 1:1
Row: 2 -> 2:2
Row: 3 -> 3:3
Row: 4 -> 4:4
Matrix after addition is:
Row: 0 -> 0:1 1:1
Row: 1 -> 0:1 1:6 2:1
Row: 2 -> 1:1 2:11 3:1
Row: 3 -> 2:1 3:16 4:1
Row: 4 -> 3:1 4:21
Add elements on CUDA device:
Initial matrix is:
Row: 0 ->
Row: 1 -> 1:1
Row: 2 -> 2:2
Row: 3 -> 3:3
Row: 4 -> 4:4
Matrix after addition is:
Row: 0 -> 0:1 1:1
Row: 1 -> 0:1 1:6 2:1
Row: 2 -> 1:1 2:11 3:1
Row: 3 -> 2:1 3:16 4:1
Row: 4 -> 3:1 4:21

◆ forAllElements() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::forAllElements ( Function &  function)

This method calls forElements for all matrix rows.

See MultidiagonalMatrix::forElements.

Template Parameters
Functionis a type of lambda function that will operate on matrix elements.
Parameters
functionis an instance of the lambda function to be called in each row.
Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forAllElementsExample()
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 0 / 1 . . . . \ -> { 0, 0, 1 }
* 0 | 2 1 . . . | -> { 0, 2, 1 }
* | 3 2 1 . . | -> { 3, 2, 1 }
* | . 3 2 1 . | -> { 3, 2, 1 }
* \ . . 3 2 1 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -2, -1, 0 }.
*/
5, // number of matrix rows
5, // number of matrix columns
{ -2, -1, 0 } ); // matrix diagonals offsets
auto f = [=] __cuda_callable__ ( int rowIdx, int localIdx, int columnIdx, double& value ) {
/***
* 'forElements' method iterates only over matrix elements lying on given subdiagonals
* and so we do not need to check anything. The element value can be expressed
* by the 'localIdx' variable, see the following figure:
*
* 0 1 2 <- localIdx values
* -------
* 0 0 / 1 . . . . \ -> { 0, 0, 1 }
* 0 | 2 1 . . . | -> { 0, 2, 1 }
* | 3 2 1 . . | -> { 3, 2, 1 }
* | . 3 2 1 . | -> { 3, 2, 1 }
* \ . . 3 2 1 / -> { 3, 2, 1 }
*
*/
value = 3 - localIdx;
};
matrix.forAllElements( f );
std::cout << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forAllElementsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forAllElementsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1
Creating matrix on CUDA device:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1

◆ forAllElements() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::forAllElements ( Function &  function) const

This method calls forElements for all matrix rows (for constant instances).

See MultidiagonalMatrix::forElements.

Template Parameters
Functionis a type of lambda function that will operate on matrix elements.
Parameters
functionis an instance of the lambda function to be called in each row.
Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forAllElementsExample()
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 0 / 1 . . . . \ -> { 0, 0, 1 }
* 0 | 2 1 . . . | -> { 0, 2, 1 }
* | 3 2 1 . . | -> { 3, 2, 1 }
* | . 3 2 1 . | -> { 3, 2, 1 }
* \ . . 3 2 1 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -2, -1, 0 }.
*/
5, // number of matrix rows
5, // number of matrix columns
{ -2, -1, 0 } ); // matrix diagonals offsets
auto f = [=] __cuda_callable__ ( int rowIdx, int localIdx, int columnIdx, double& value ) {
/***
* 'forElements' method iterates only over matrix elements lying on given subdiagonals
* and so we do not need to check anything. The element value can be expressed
* by the 'localIdx' variable, see the following figure:
*
* 0 1 2 <- localIdx values
* -------
* 0 0 / 1 . . . . \ -> { 0, 0, 1 }
* 0 | 2 1 . . . | -> { 0, 2, 1 }
* | 3 2 1 . . | -> { 3, 2, 1 }
* | . 3 2 1 . | -> { 3, 2, 1 }
* \ . . 3 2 1 / -> { 3, 2, 1 }
*
*/
value = 3 - localIdx;
};
matrix.forAllElements( f );
std::cout << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forAllElementsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forAllElementsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1
Creating matrix on CUDA device:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1

◆ forAllRows() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::forAllRows ( Function &&  function)

Method for parallel iteration over all matrix rows.

In each row, given lambda function is performed. Each row is processed by at most one thread unlike the method MultidiagonalMatrix::forAllElements where more than one thread can be mapped to each row.

Template Parameters
Functionis type of the lambda function.
Parameters
functionis an instance of the lambda function to be called for each row.
auto function = [] __cuda_callable__ ( RowView& row ) mutable { ... };
typename ViewType::RowView RowView
Type for accessing matrix rows.
Definition: MultidiagonalMatrix.h:136

RowView represents matrix row - see TNL::Matrices::MultidiagonalMatrix::RowView.

Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forRowsExample()
{
using RowView = typename MatrixType::RowView;
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 2 -1. . . . \ -> { 0, 0, 1 }
* | -1 2 -1 . . | -> { 0, 2, 1 }
* | . -1 2 -1. . | -> { 3, 2, 1 }
* | . . -1 2 -1 | -> { 3, 2, 1 }
* \ . . . -1. 2 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -1, 0, 1 }.
*/
const int size = 5;
MatrixType matrix(
size, // number of matrix rows
size, // number of matrix columns
{ -1, -0, 1 } ); // matrix diagonals offsets
auto f = [=] __cuda_callable__ ( RowView& row ) {
const int& rowIdx = row.getRowIndex();
if( rowIdx > 0 )
row.setElement( 0, -1.0 ); // elements below the diagonal
row.setElement( 1, 2.0 ); // elements on the diagonal
if( rowIdx < size - 1 ) // elements above the diagonal
row.setElement( 2, -1.0 );
};
matrix.forAllRows( f );
std::cout << matrix << std::endl;
/***
* Compute sum of elements in each row and store it into a vector.
*/
auto sum_view = sum_vector.getView();
matrix.forAllRows( [=] __cuda_callable__ ( RowView& row ) mutable {
double sum( 0.0 );
for( auto element : row )
sum += TNL::abs( element.value() );
sum_view[ row.getRowIndex() ] = sum;
} );
std::cout << "Sums in matrix rows = " << sum_vector << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forRowsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forRowsExample< TNL::Devices::Cuda >();
#endif
}
__cuda_callable__ T abs(const T &n)
This function returns absolute value of given number n.
Definition: Math.h:91
Structure for specifying type of sparse matrix.
Definition: MatrixType.h:21
Output
Creating matrix on host:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Sums in matrix rows = [ 3, 4, 4, 4, 3 ]
Creating matrix on CUDA device:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Sums in matrix rows = [ 3, 4, 4, 4, 3 ]

◆ forAllRows() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::forAllRows ( Function &&  function) const

Method for parallel iteration over all matrix rows for constant instances.

In each row, given lambda function is performed. Each row is processed by at most one thread unlike the method MultidiagonalMatrix::forAllElements where more than one thread can be mapped to each row.

Template Parameters
Functionis type of the lambda function.
Parameters
functionis an instance of the lambda function to be called for each row.
auto function = [] __cuda_callable__ ( RowView& row ) { ... };

RowView represents matrix row - see TNL::Matrices::MultidiagonalMatrix::RowView.

Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forRowsExample()
{
using RowView = typename MatrixType::RowView;
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 2 -1. . . . \ -> { 0, 0, 1 }
* | -1 2 -1 . . | -> { 0, 2, 1 }
* | . -1 2 -1. . | -> { 3, 2, 1 }
* | . . -1 2 -1 | -> { 3, 2, 1 }
* \ . . . -1. 2 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -1, 0, 1 }.
*/
const int size = 5;
MatrixType matrix(
size, // number of matrix rows
size, // number of matrix columns
{ -1, -0, 1 } ); // matrix diagonals offsets
auto f = [=] __cuda_callable__ ( RowView& row ) {
const int& rowIdx = row.getRowIndex();
if( rowIdx > 0 )
row.setElement( 0, -1.0 ); // elements below the diagonal
row.setElement( 1, 2.0 ); // elements on the diagonal
if( rowIdx < size - 1 ) // elements above the diagonal
row.setElement( 2, -1.0 );
};
matrix.forAllRows( f );
std::cout << matrix << std::endl;
/***
* Compute sum of elements in each row and store it into a vector.
*/
auto sum_view = sum_vector.getView();
matrix.forAllRows( [=] __cuda_callable__ ( RowView& row ) mutable {
double sum( 0.0 );
for( auto element : row )
sum += TNL::abs( element.value() );
sum_view[ row.getRowIndex() ] = sum;
} );
std::cout << "Sums in matrix rows = " << sum_vector << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forRowsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forRowsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Sums in matrix rows = [ 3, 4, 4, 4, 3 ]
Creating matrix on CUDA device:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Sums in matrix rows = [ 3, 4, 4, 4, 3 ]

◆ forElements() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::forElements ( IndexType  begin,
IndexType  end,
Function &  function 
)

Method for iteration over matrix rows for non-constant instances.

Template Parameters
Functionis type of lambda function that will operate on matrix elements. It is should have form like
auto function = [=] __cuda_callable__ ( IndexType rowIdx, IndexType localIdx, IndexType columnIdx, const RealType& value )
{ ... };
Index IndexType
The type used for matrix elements indexing.
Definition: MultidiagonalMatrix.h:96
Definition: Real.h:17

where

rowIdx is an index of the matrix row.

localIdx parameter is a rank of the non-zero element in given row. It is also, in fact, index of the matrix subdiagonal.

columnIdx is a column index of the matrix element.

value is a reference to the matrix element value. It can be used even for changing the matrix element value.

Parameters
begindefines beginning of the range [begin,end) of rows to be processed.
enddefines ending of the range [begin,end) of rows to be processed.
functionis an instance of the lambda function to be called in each row.
Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forRowsExample()
{
using RowView = typename MatrixType::RowView;
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 2 -1. . . . \ -> { 0, 0, 1 }
* | -1 2 -1 . . | -> { 0, 2, 1 }
* | . -1 2 -1. . | -> { 3, 2, 1 }
* | . . -1 2 -1 | -> { 3, 2, 1 }
* \ . . . -1. 2 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -1, 0, 1 }.
*/
const int size = 5;
MatrixType matrix(
size, // number of matrix rows
size, // number of matrix columns
{ -1, -0, 1 } ); // matrix diagonals offsets
auto f = [=] __cuda_callable__ ( RowView& row ) {
const int& rowIdx = row.getRowIndex();
if( rowIdx > 0 )
row.setElement( 0, -1.0 ); // elements below the diagonal
row.setElement( 1, 2.0 ); // elements on the diagonal
if( rowIdx < size - 1 ) // elements above the diagonal
row.setElement( 2, -1.0 );
};
matrix.forAllRows( f );
std::cout << matrix << std::endl;
/***
* Compute sum of elements in each row and store it into a vector.
*/
auto sum_view = sum_vector.getView();
matrix.forAllRows( [=] __cuda_callable__ ( RowView& row ) mutable {
double sum( 0.0 );
for( auto element : row )
sum += TNL::abs( element.value() );
sum_view[ row.getRowIndex() ] = sum;
} );
std::cout << "Sums in matrix rows = " << sum_vector << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forRowsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forRowsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Sums in matrix rows = [ 3, 4, 4, 4, 3 ]
Creating matrix on CUDA device:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Sums in matrix rows = [ 3, 4, 4, 4, 3 ]

◆ forElements() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::forElements ( IndexType  begin,
IndexType  end,
Function &  function 
) const

Method for iteration over matrix rows for constant instances.

Template Parameters
Functionis type of lambda function that will operate on matrix elements. It is should have form like
auto function = [=] __cuda_callble__ ( IndexType rowIdx, IndexType localIdx, IndexType columnIdx, const RealType& value )
{ ... };

where

rowIdx is an index of the matrix row.

localIdx parameter is a rank of the non-zero element in given row. It is also, in fact, index of the matrix subdiagonal.

columnIdx is a column index of the matrix element.

value is the matrix element value.

Parameters
begindefines beginning of the range [begin,end) of rows to be processed.
enddefines ending of the range [begin,end) of rows to be processed.
functionis an instance of the lambda function to be called in each row.
Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forRowsExample()
{
using RowView = typename MatrixType::RowView;
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 2 -1. . . . \ -> { 0, 0, 1 }
* | -1 2 -1 . . | -> { 0, 2, 1 }
* | . -1 2 -1. . | -> { 3, 2, 1 }
* | . . -1 2 -1 | -> { 3, 2, 1 }
* \ . . . -1. 2 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -1, 0, 1 }.
*/
const int size = 5;
MatrixType matrix(
size, // number of matrix rows
size, // number of matrix columns
{ -1, -0, 1 } ); // matrix diagonals offsets
auto f = [=] __cuda_callable__ ( RowView& row ) {
const int& rowIdx = row.getRowIndex();
if( rowIdx > 0 )
row.setElement( 0, -1.0 ); // elements below the diagonal
row.setElement( 1, 2.0 ); // elements on the diagonal
if( rowIdx < size - 1 ) // elements above the diagonal
row.setElement( 2, -1.0 );
};
matrix.forAllRows( f );
std::cout << matrix << std::endl;
/***
* Compute sum of elements in each row and store it into a vector.
*/
auto sum_view = sum_vector.getView();
matrix.forAllRows( [=] __cuda_callable__ ( RowView& row ) mutable {
double sum( 0.0 );
for( auto element : row )
sum += TNL::abs( element.value() );
sum_view[ row.getRowIndex() ] = sum;
} );
std::cout << "Sums in matrix rows = " << sum_vector << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forRowsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forRowsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Sums in matrix rows = [ 3, 4, 4, 4, 3 ]
Creating matrix on CUDA device:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Sums in matrix rows = [ 3, 4, 4, 4, 3 ]

◆ forRows() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::forRows ( IndexType  begin,
IndexType  end,
Function &&  function 
)

Method for parallel iteration over matrix rows from interval [ begin, end).

In each row, given lambda function is performed. Each row is processed by at most one thread unlike the method MultidiagonalMatrix::forElements where more than one thread can be mapped to each row.

Template Parameters
Functionis type of the lambda function.
Parameters
begindefines beginning of the range [ begin,end ) of rows to be processed.
enddefines ending of the range [ begin, end ) of rows to be processed.
functionis an instance of the lambda function to be called for each row.
auto function = [] __cuda_callable__ ( RowView& row ) mutable { ... };

RowView represents matrix row - see TNL::Matrices::MultidiagonalMatrix::RowView.

Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forRowsExample()
{
using RowView = typename MatrixType::RowView;
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 2 -1. . . . \ -> { 0, 0, 1 }
* | -1 2 -1 . . | -> { 0, 2, 1 }
* | . -1 2 -1. . | -> { 3, 2, 1 }
* | . . -1 2 -1 | -> { 3, 2, 1 }
* \ . . . -1. 2 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -1, 0, 1 }.
*/
const int size = 5;
MatrixType matrix(
size, // number of matrix rows
size, // number of matrix columns
{ -1, -0, 1 } ); // matrix diagonals offsets
auto f = [=] __cuda_callable__ ( RowView& row ) {
const int& rowIdx = row.getRowIndex();
if( rowIdx > 0 )
row.setElement( 0, -1.0 ); // elements below the diagonal
row.setElement( 1, 2.0 ); // elements on the diagonal
if( rowIdx < size - 1 ) // elements above the diagonal
row.setElement( 2, -1.0 );
};
matrix.forAllRows( f );
std::cout << matrix << std::endl;
/***
* Compute sum of elements in each row and store it into a vector.
*/
auto sum_view = sum_vector.getView();
matrix.forAllRows( [=] __cuda_callable__ ( RowView& row ) mutable {
double sum( 0.0 );
for( auto element : row )
sum += TNL::abs( element.value() );
sum_view[ row.getRowIndex() ] = sum;
} );
std::cout << "Sums in matrix rows = " << sum_vector << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forRowsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forRowsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Sums in matrix rows = [ 3, 4, 4, 4, 3 ]
Creating matrix on CUDA device:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Sums in matrix rows = [ 3, 4, 4, 4, 3 ]

◆ forRows() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::forRows ( IndexType  begin,
IndexType  end,
Function &&  function 
) const

Method for parallel iteration over matrix rows from interval [ begin, end) for constant instances.

In each row, given lambda function is performed. Each row is processed by at most one thread unlike the method MultidiagonalMatrix::forElements where more than one thread can be mapped to each row.

Template Parameters
Functionis type of the lambda function.
Parameters
begindefines beginning of the range [ begin,end ) of rows to be processed.
enddefines ending of the range [ begin, end ) of rows to be processed.
functionis an instance of the lambda function to be called for each row.
auto function = [] __cuda_callable__ ( RowView& row ) { ... };

RowView represents matrix row - see TNL::Matrices::MultidiagonalMatrix::RowView.

Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forRowsExample()
{
using RowView = typename MatrixType::RowView;
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 2 -1. . . . \ -> { 0, 0, 1 }
* | -1 2 -1 . . | -> { 0, 2, 1 }
* | . -1 2 -1. . | -> { 3, 2, 1 }
* | . . -1 2 -1 | -> { 3, 2, 1 }
* \ . . . -1. 2 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -1, 0, 1 }.
*/
const int size = 5;
MatrixType matrix(
size, // number of matrix rows
size, // number of matrix columns
{ -1, -0, 1 } ); // matrix diagonals offsets
auto f = [=] __cuda_callable__ ( RowView& row ) {
const int& rowIdx = row.getRowIndex();
if( rowIdx > 0 )
row.setElement( 0, -1.0 ); // elements below the diagonal
row.setElement( 1, 2.0 ); // elements on the diagonal
if( rowIdx < size - 1 ) // elements above the diagonal
row.setElement( 2, -1.0 );
};
matrix.forAllRows( f );
std::cout << matrix << std::endl;
/***
* Compute sum of elements in each row and store it into a vector.
*/
auto sum_view = sum_vector.getView();
matrix.forAllRows( [=] __cuda_callable__ ( RowView& row ) mutable {
double sum( 0.0 );
for( auto element : row )
sum += TNL::abs( element.value() );
sum_view[ row.getRowIndex() ] = sum;
} );
std::cout << "Sums in matrix rows = " << sum_vector << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forRowsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forRowsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Sums in matrix rows = [ 3, 4, 4, 4, 3 ]
Creating matrix on CUDA device:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Sums in matrix rows = [ 3, 4, 4, 4, 3 ]

◆ getCompressedRowLengths()

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 >::getCompressedRowLengths ( Vector &  rowLengths) const

Computes number of non-zeros in each row.

Parameters
rowLengthsis a vector into which the number of non-zeros in each row will be stored.
Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void laplaceOperatorMatrix()
{
const int gridSize( 4 );
const int matrixSize = gridSize * gridSize;
matrixSize, // number of rows
matrixSize, // number of columns
{ - gridSize, -1, 0, 1, gridSize } // diagonals offsets
);
matrix.setElements( {
{ 0.0, 0.0, 1.0 }, // set matrix elements corresponding to boundary grid nodes
{ 0.0, 0.0, 1.0 }, // and Dirichlet boundary conditions, i.e. 1 on the main diagonal
{ 0.0, 0.0, 1.0 }, // which is the third one
{ 0.0, 0.0, 1.0 },
{ 0.0, 0.0, 1.0 },
{ -1.0, -1.0, 4.0, -1.0, -1.0 }, // set matrix elements corresponding to inner grid nodes, i.e. 4 on the main diagonal
{ -1.0, -1.0, 4.0, -1.0, -1.0 }, // (the third one) and -1 to the other sub-diagonals
{ 0.0, 0.0, 1.0 },
{ 0.0, 0.0, 1.0 },
{ -1.0, -1.0, 4.0, -1.0, -1.0 },
{ -1.0, -1.0, 4.0, -1.0, -1.0 },
{ 0.0, 0.0, 1.0 },
{ 0.0, 0.0, 1.0 },
{ 0.0, 0.0, 1.0 },
{ 0.0, 0.0, 1.0 },
{ 0.0, 0.0, 1.0 }
} );
matrix.getCompressedRowLengths( rowLengths );
std::cout << "Laplace operator matrix: " << std::endl << matrix << std::endl;
std::cout << "Compressed row lengths: " << rowLengths << 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
}
void setElements(const std::initializer_list< std::initializer_list< ListReal > > &data)
Set matrix elements from an initializer list.
Definition: MultidiagonalMatrix.hpp:230
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
Compressed row lengths: [ 1, 1, 1, 1, 1, 5, 5, 1, 1, 5, 5, 1, 1, 1, 1, 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
Compressed row lengths: [ 1, 1, 1, 1, 1, 5, 5, 1, 1, 5, 5, 1, 1, 1, 1, 1 ]

◆ 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

Returns a non-modifiable view of the multidiagonal matrix.

See MultidiagonalMatrixView.

Returns
multidiagonal matrix view.

◆ getDiagonalsCount()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
Index TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::getDiagonalsCount

Returns number of diagonals.

Returns
Number of diagonals.

◆ getDiagonalsOffsets()

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

Returns vector with diagonals offsets.

Returns
vector with diagonals offsets.

◆ getElement()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
__cuda_callable__ Real TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::getElement ( IndexType  row,
IndexType  column 
) const

Returns value of matrix element at position given by its row and column index.

This method can be called from the host system (CPU) no matter where the matrix is allocated. If the matrix is allocated on GPU this method can be called even from device kernels. If the matrix is allocated in GPU device this method is called from CPU, it transfers values of each matrix element separately and so the performance is very low. For higher performance see. MultidiagonalMatrix::getRow or MultidiagonalMatrix::forElements and MultidiagonalMatrix::forAllElements.

Parameters
rowis a row index of the matrix element.
columni a column index of the matrix element.
Returns
value of given matrix element.
Example
#include <iostream>
#include <iomanip>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void getElements()
{
const int matrixSize( 5 );
matrixSize, // number of matrix columns
{ -1, 0, 1 }, // matrix diagonals offsets
{ // matrix elements definition
{ 0.0, 2.0, -1.0 },
{ -1.0, 2.0, -1.0 },
{ -1.0, 2.0, -1.0 },
{ -1.0, 2.0, -1.0 },
{ -1.0, 2.0, 0.0 }
} );
for( int i = 0; i < matrixSize; i++ )
{
for( int j = 0; j < matrixSize; j++ )
std::cout << std::setw( 5 ) << matrix.getElement( i, j );
}
}
int main( int argc, char* argv[] )
{
std::cout << "Get elements on host:" << std::endl;
getElements< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Get elements on CUDA device:" << std::endl;
getElements< TNL::Devices::Cuda >();
#endif
}
T setw(T... args)
Output
Get elements on host:
2 -1 0 0 0
-1 2 -1 0 0
0 -1 2 -1 0
0 0 -1 2 -1
0 0 0 -1 2
Get elements on CUDA device:
2 -1 0 0 0
-1 2 -1 0 0
0 -1 2 -1 0
0 0 -1 2 -1
0 0 0 -1 2

◆ getIndexer() [1/2]

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

This method returns matrix elements indexer used by this matrix.

Returns
non-constant reference to the indexer.

◆ getIndexer() [2/2]

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

This method returns matrix elements indexer used by this matrix.

Returns
constant reference to the indexer.

◆ getNonzeroElementsCount()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
Index TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::getNonzeroElementsCount
overridevirtual

Returns number of non-zero matrix elements.

This method really counts the non-zero matrix elements and so it returns zero for matrix having all allocated elements set to zero.

Returns
number of non-zero matrix elements.

Reimplemented from TNL::Matrices::Matrix< Real, Device, Index, RealAllocator >.

◆ getPaddingIndex()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
__cuda_callable__ Index TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::getPaddingIndex

Returns padding index denoting padding zero elements.

These elements are used for efficient data alignment in memory.

Returns
value of the padding index.

◆ getRow() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
__cuda_callable__ auto TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::getRow ( IndexType  rowIdx)

Non-constant getter of simple structure for accessing given matrix row.

Parameters
rowIdxis matrix row index.
Returns
RowView for accessing given matrix row.
Example
#include <iostream>
#include <TNL/Algorithms/parallelFor.h>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
#include <TNL/Pointers/SharedPointer.h>
template< typename Device >
void getRowExample()
{
const int matrixSize( 5 );
auto diagonalsOffsets = { -1, 0, 1 }; // Variadic templates in SharedPointer
// constructor do not recognize initializer
// list so we give it a hint.
matrixSize, // number of matrix rows
matrixSize, // number of matrix columns
diagonalsOffsets );
MatrixType* matrix_device = &matrix.template modifyData< Device >();
auto f = [=] __cuda_callable__ ( int rowIdx ) mutable {
auto row = matrix_device->getRow( rowIdx );
if( rowIdx > 0 )
row.setElement( 0, -1.0 ); // elements below the diagonal
row.setElement( 1, 2.0 ); // elements on the diagonal
if( rowIdx < matrixSize - 1 ) // elements above the diagonal
row.setElement( 2, -1.0 );
};
/***
* For the case when Device is CUDA device we need to synchronize smart
* pointers. To avoid this you may use MultidiagonalMatrixView. See
* MultidiagonalMatrixView::getRow example for details.
*/
TNL::Pointers::synchronizeSmartPointersOnDevice< Device >();
/***
* Set the matrix elements.
*/
TNL::Algorithms::parallelFor< Device >( 0, matrix->getRows(), f );
std::cout << std::endl << *matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Getting matrix rows on host: " << std::endl;
getRowExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Getting matrix rows on CUDA device: " << std::endl;
getRowExample< TNL::Devices::Cuda >();
#endif
}
Cross-device shared smart pointer.
Definition: SharedPointer.h:49
Output
Getting matrix rows on host:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2
Getting matrix rows on CUDA device:
Row: 0 -> 0:2 1:-1
Row: 1 -> 0:-1 1:2 2:-1
Row: 2 -> 1:-1 2:2 3:-1
Row: 3 -> 2:-1 3:2 4:-1
Row: 4 -> 3:-1 4:2

See MultidiagonalMatrixRowView.

◆ getRow() [2/2]

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

Constant getter of simple structure for accessing given matrix row.

Parameters
rowIdxis matrix row index.
Returns
RowView for accessing given matrix row.
Example
#include <iostream>
#include <functional>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
#include <TNL/Pointers/SharedPointer.h>
template< typename Device >
void getRowExample()
{
const int matrixSize = 5;
auto diagonalsOffsets = { -2, -1, 0 };
matrixSize, // number of matrix rows
matrixSize, // number of matrix columns
diagonalsOffsets );
matrix->setElements(
{ { 0.0, 0.0, 1.0 },
{ 0.0, 2.0, 1.0 },
{ 3.0, 2.0, 1.0 },
{ 3.0, 2.0, 1.0 },
{ 3.0, 2.0, 1.0 } } );
/***
* Fetch lambda function returns diagonal element in each row.
*/
const MatrixType* matrix_device = &matrix.template getData< Device >();
auto fetch = [=] __cuda_callable__ ( int rowIdx ) mutable -> double {
auto row = matrix_device->getRow( rowIdx );
return row.getValue( 2 ); // get value from subdiagonal with index 2, i.e. the main diagonal
};
/***
* For the case when Device is CUDA device we need to synchronize smart
* pointers. To avoid this you may use MultidiagonalMatrixView. See
* MultidiagonalMatrixView::getConstRow example for details.
*/
TNL::Pointers::synchronizeSmartPointersOnDevice< Device >();
/***
* Compute the matrix trace.
*/
int trace = TNL::Algorithms::reduce< Device >( 0, matrix->getRows(), fetch, std::plus<>{}, 0 );
std::cout << "Matrix reads as: " << std::endl << *matrix << std::endl;
std::cout << "Matrix trace is: " << trace << "." << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Getting matrix rows on host: " << std::endl;
getRowExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Getting matrix rows on CUDA device: " << std::endl;
getRowExample< TNL::Devices::Cuda >();
#endif
}
Output
Getting matrix rows on host:
Matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1
Matrix trace is: 5.
Getting matrix rows on CUDA device:
Matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1
Matrix trace is: 5.

See MultidiagonalMatrixRowView.

◆ getRowCapacities()

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 >::getRowCapacities ( Vector &  rowCapacities) const

Compute capacities of all rows.

The row capacities are not stored explicitly and must be computed.

Parameters
rowCapacitiesis a vector where the row capacities will be stored.

◆ getSerializationType()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
std::string TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::getSerializationType
static

Returns string with serialization type.

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

Returns
String with the serialization type.
Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void getSerializationTypeExample()
{
std::cout << "Matrix type is: " << matrix.getSerializationType();
}
int main( int argc, char* argv[] )
{
std::cout << "Get serialization type on CPU ... " << std::endl;
getSerializationTypeExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Get serialization type on CUDA GPU ... " << std::endl;
getSerializationTypeExample< TNL::Devices::Cuda >();
#endif
}
static std::string getSerializationType()
Returns string with serialization type.
Definition: MultidiagonalMatrix.hpp:107
Output
Get serialization type on CPU ...
Matrix type is: Matrices::MultidiagonalMatrix< double, [any_device], int, RowMajorOrder, [any_allocator], [any_allocator] >Get serialization type on CUDA GPU ...
Matrix type is: Matrices::MultidiagonalMatrix< double, [any_device], int, ColumnMajorOrder, [any_allocator], [any_allocator] >

◆ getSerializationTypeVirtual()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
std::string TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::getSerializationTypeVirtual
overridevirtual

Returns string with serialization type.

See MultidiagonalMatrix::getSerializationType.

Returns
String with the serialization type.
Example
#include <iostream>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void getSerializationTypeExample()
{
std::cout << "Matrix type is: " << matrix.getSerializationType();
}
int main( int argc, char* argv[] )
{
std::cout << "Get serialization type on CPU ... " << std::endl;
getSerializationTypeExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Get serialization type on CUDA GPU ... " << std::endl;
getSerializationTypeExample< TNL::Devices::Cuda >();
#endif
}
Output
Get serialization type on CPU ...
Matrix type is: Matrices::MultidiagonalMatrix< double, [any_device], int, RowMajorOrder, [any_allocator], [any_allocator] >Get serialization type on CUDA GPU ...
Matrix type is: Matrices::MultidiagonalMatrix< double, [any_device], int, ColumnMajorOrder, [any_allocator], [any_allocator] >

Reimplemented from TNL::Object.

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

◆ isSymmetric()

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 >>
static constexpr bool TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::isSymmetric ( )
inlinestaticconstexpr

This is only for compatibility with sparse matrices.

Returns
false.

◆ 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::Matrices::Matrix< Real, Device, Index, RealAllocator >.

◆ operator!=()

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_ >
bool TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::operator!= ( const MultidiagonalMatrix< Real_, Device_, Index_, Organization_, RealAllocator_, IndexAllocator_ > &  matrix) const

Comparison operator with 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.
Returns
true if both matrices are NOT identical and false otherwise.

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

◆ operator==()

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_ >
bool TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::operator== ( const MultidiagonalMatrix< Real_, Device_, Index_, Organization_, RealAllocator_, IndexAllocator_ > &  matrix) const

Comparison operator with 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.
Returns
true if both matrices are identical and false otherwise.

◆ print()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::print ( std::ostream str) const
overridevirtual

Method for printing the matrix to output stream.

Parameters
stris the output stream.

Reimplemented from TNL::Matrices::Matrix< Real, Device, Index, RealAllocator >.

◆ reduceAllRows() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::reduceAllRows ( Fetch &  fetch,
Reduce &  reduce,
Keep &  keep,
const FetchReal &  identity 
)

Method for performing general reduction on all matrix rows.

Template Parameters
Fetchis a type of lambda function for data fetch declared as
auto fetch = [=] __cuda_callable__ ( IndexType rowIdx, IndexType& columnIdx, RealType& elementValue ) -> FetchValue { ...
};

The return type of this lambda can be any non void.

Template Parameters
Reduceis a type of lambda function for reduction declared as
auto reduce = [=] __cuda_callable__ ( const FetchValue& v1, const FetchValue& v2 ) -> FetchValue { ... };
Template Parameters
Keepis a type of lambda function for storing results of reduction in each row. It is declared as
auto keep = [=] __cuda_callable__ ( IndexType rowIdx, const RealType& value ) { ... };
Template Parameters
FetchValueis type returned by the Fetch lambda function.
Parameters
fetchis an instance of lambda function for data fetch.
reduceis an instance of lambda function for reduction.
keepin an instance of lambda function for storing results.
identityis the identity element for the reduction operation, i.e. element which does not change the result of the reduction.
Example
#include <iostream>
#include <iomanip>
#include <functional>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 0 / 1 . . . . \ -> { 0, 0, 1 }
* 0 | 2 1 . . . | -> { 0, 2, 1 }
* | 3 2 1 . . | -> { 3, 2, 1 }
* | . 3 2 1 . | -> { 3, 2, 1 }
* \ . . 3 2 1 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -2, -1, 0 }.
*/
5, // number of matrix columns
{ -2, -1, 0 }, // diagonals offsets
{ { 0, 0, 1 }, // matrix elements
{ 0, 2, 1 },
{ 3, 2, 1 },
{ 3, 2, 1 },
{ 3, 2, 1 } } );
/***
* Find largest element in each row.
*/
TNL::Containers::Vector< double, Device > rowMax( matrix.getRows() );
/***
* Prepare vector view for lambdas.
*/
auto rowMaxView = rowMax.getView();
/***
* Fetch lambda just returns absolute value of matrix elements.
*/
auto fetch = [=] __cuda_callable__ ( int rowIdx, int columnIdx, const double& value ) -> double {
return TNL::abs( value );
};
/***
* Reduce lambda return maximum of given values.
*/
auto reduce = [=] __cuda_callable__ ( double& a, const double& b ) -> double {
return TNL::max( a, b );
};
/***
* Keep lambda store the largest value in each row to the vector rowMax.
*/
auto keep = [=] __cuda_callable__ ( int rowIdx, const double& value ) mutable {
rowMaxView[ rowIdx ] = value;
};
/***
* Compute the largest values in each row.
*/
matrix.reduceAllRows( fetch, reduce, keep, std::numeric_limits< double >::lowest() );
std::cout << "The matrix reads as: " << std::endl << matrix << std::endl;
std::cout << "Max. elements in rows are: " << rowMax << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Rows reduction on host:" << std::endl;
reduceAllRows< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Rows reduction on CUDA device:" << std::endl;
reduceAllRows< TNL::Devices::Cuda >();
#endif
}
void reduceAllRows(Fetch &fetch, Reduce &reduce, Keep &keep, const FetchReal &identity)
Method for performing general reduction on all matrix rows.
Definition: MultidiagonalMatrix.hpp:536
Result reduce(Index begin, Index end, Fetch &&fetch, Reduction &&reduction, const Result &identity)
reduce implements (parallel) reduction for vectors and arrays.
Definition: reduce.h:71
constexpr ResultType max(const T1 &a, const T2 &b)
This function returns maximum of two numbers.
Definition: Math.h:65
Output
Rows reduction on host:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1
Max. elements in rows are: [ 1, 2, 3, 3, 3 ]
Rows reduction on CUDA device:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1
Max. elements in rows are: [ 1, 2, 3, 3, 3 ]

◆ reduceAllRows() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::reduceAllRows ( Fetch &  fetch,
Reduce &  reduce,
Keep &  keep,
const FetchReal &  identity 
) const

Method for performing general reduction on all matrix rows for constant instances.

Template Parameters
Fetchis a type of lambda function for data fetch declared as
auto fetch = [=] __cuda_callable__ ( IndexType rowIdx, IndexType& columnIdx, RealType& elementValue ) -> FetchValue { ...
};

The return type of this lambda can be any non void.

Template Parameters
Reduceis a type of lambda function for reduction declared as
auto reduce = [=] __cuda_callable__ ( const FetchValue& v1, const FetchValue& v2 ) -> FetchValue { ... };
Template Parameters
Keepis a type of lambda function for storing results of reduction in each row. It is declared as
auto keep = [=] __cuda_callable__ ( IndexType rowIdx, const RealType& value ) { ... };
Template Parameters
FetchValueis type returned by the Fetch lambda function.
Parameters
fetchis an instance of lambda function for data fetch.
reduceis an instance of lambda function for reduction.
keepin an instance of lambda function for storing results.
identityis the identity element for the reduction operation, i.e. element which does not change the result of the reduction.
Example
#include <iostream>
#include <iomanip>
#include <functional>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 0 / 1 . . . . \ -> { 0, 0, 1 }
* 0 | 2 1 . . . | -> { 0, 2, 1 }
* | 3 2 1 . . | -> { 3, 2, 1 }
* | . 3 2 1 . | -> { 3, 2, 1 }
* \ . . 3 2 1 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -2, -1, 0 }.
*/
5, // number of matrix columns
{ -2, -1, 0 }, // diagonals offsets
{ { 0, 0, 1 }, // matrix elements
{ 0, 2, 1 },
{ 3, 2, 1 },
{ 3, 2, 1 },
{ 3, 2, 1 } } );
/***
* Find largest element in each row.
*/
TNL::Containers::Vector< double, Device > rowMax( matrix.getRows() );
/***
* Prepare vector view for lambdas.
*/
auto rowMaxView = rowMax.getView();
/***
* Fetch lambda just returns absolute value of matrix elements.
*/
auto fetch = [=] __cuda_callable__ ( int rowIdx, int columnIdx, const double& value ) -> double {
return TNL::abs( value );
};
/***
* Reduce lambda return maximum of given values.
*/
auto reduce = [=] __cuda_callable__ ( double& a, const double& b ) -> double {
return TNL::max( a, b );
};
/***
* Keep lambda store the largest value in each row to the vector rowMax.
*/
auto keep = [=] __cuda_callable__ ( int rowIdx, const double& value ) mutable {
rowMaxView[ rowIdx ] = value;
};
/***
* Compute the largest values in each row.
*/
matrix.reduceAllRows( fetch, reduce, keep, std::numeric_limits< double >::lowest() );
std::cout << "The matrix reads as: " << std::endl << matrix << std::endl;
std::cout << "Max. elements in rows are: " << rowMax << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Rows reduction on host:" << std::endl;
reduceAllRows< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Rows reduction on CUDA device:" << std::endl;
reduceAllRows< TNL::Devices::Cuda >();
#endif
}
Output
Rows reduction on host:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1
Max. elements in rows are: [ 1, 2, 3, 3, 3 ]
Rows reduction on CUDA device:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1
Max. elements in rows are: [ 1, 2, 3, 3, 3 ]

◆ reduceRows() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::reduceRows ( IndexType  begin,
IndexType  end,
Fetch &  fetch,
Reduce &  reduce,
Keep &  keep,
const FetchReal &  identity 
)

Method for performing general reduction on matrix rows.

Template Parameters
Fetchis a type of lambda function for data fetch declared as
auto fetch = [=] __cuda_callable__ ( IndexType rowIdx, IndexType& columnIdx, RealType& elementValue ) -> FetchValue { ...
};

The return type of this lambda can be any non void.

Template Parameters
Reduceis a type of lambda function for reduction declared as
auto reduce = [=] __cuda_callable__ ( const FetchValue& v1, const FetchValue& v2 ) -> FetchValue { ... };
Template Parameters
Keepis a type of lambda function for storing results of reduction in each row. It is declared as
auto keep = [=] __cuda_callable__ ( IndexType rowIdx, const RealType& value ) { ... };
Template Parameters
FetchValueis type returned by the Fetch lambda function.
Parameters
begindefines beginning of the range [begin,end) of rows to be processed.
enddefines ending of the range [begin,end) of rows to be processed.
fetchis an instance of lambda function for data fetch.
reduceis an instance of lambda function for reduction.
keepin an instance of lambda function for storing results.
identityis the identity element for the reduction operation, i.e. element which does not change the result of the reduction.
Example
#include <iostream>
#include <iomanip>
#include <functional>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void reduceRows()
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 0 / 1 . . . . \ -> { 0, 0, 1 }
* 0 | 2 1 . . . | -> { 0, 2, 1 }
* | 3 2 1 . . | -> { 3, 2, 1 }
* | . 3 2 1 . | -> { 3, 2, 1 }
* \ . . 3 2 1 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -2, -1, 0 }.
*/
5, // number of matrix columns
{ -2, -1, 0 }, // diagonals offsets
{ { 0, 0, 1 }, // matrix elements
{ 0, 2, 1 },
{ 3, 2, 1 },
{ 3, 2, 1 },
{ 3, 2, 1 } } );
/***
* Find largest element in each row.
*/
TNL::Containers::Vector< double, Device > rowMax( matrix.getRows() );
/***
* Prepare vector view for lambdas.
*/
auto rowMaxView = rowMax.getView();
/***
* Fetch lambda just returns absolute value of matrix elements.
*/
auto fetch = [=] __cuda_callable__ ( int rowIdx, int columnIdx, const double& value ) -> double {
return TNL::abs( value );
};
/***
* Reduce lambda return maximum of given values.
*/
auto reduce = [=] __cuda_callable__ ( const double& a, const double& b ) -> double {
return TNL::max( a, b );
};
/***
* Keep lambda store the largest value in each row to the vector rowMax.
*/
auto keep = [=] __cuda_callable__ ( int rowIdx, const double& value ) mutable {
rowMaxView[ rowIdx ] = value;
};
/***
* Compute the largest values in each row.
*/
matrix.reduceRows( 0, matrix.getRows(), fetch, reduce, keep, std::numeric_limits< double >::lowest() );
std::cout << "The matrix reads as: " << std::endl << matrix << std::endl;
std::cout << "Max. elements in rows are: " << rowMax << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Rows reduction on host:" << std::endl;
reduceRows< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Rows reduction on CUDA device:" << std::endl;
reduceRows< TNL::Devices::Cuda >();
#endif
}
void reduceRows(IndexType begin, IndexType end, Fetch &fetch, Reduce &reduce, Keep &keep, const FetchReal &identity)
Method for performing general reduction on matrix rows.
Definition: MultidiagonalMatrix.hpp:501
Output
Rows reduction on host:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1
Max. elements in rows are: [ 1, 2, 3, 3, 3 ]
Rows reduction on CUDA device:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1
Max. elements in rows are: [ 1, 2, 3, 3, 3 ]

◆ reduceRows() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::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 Parameters
Fetchis a type of lambda function for data fetch declared as
auto fetch = [=] __cuda_callable__ ( IndexType rowIdx, IndexType& columnIdx, RealType& elementValue ) -> FetchValue { ...
};

The return type of this lambda can be any non void.

Template Parameters
Reduceis a type of lambda function for reduction declared as
auto reduce = [=] __cuda_callbale__ ( const FetchValue& v1, const FetchValue& v2 ) -> FetchValue { ... };
Template Parameters
Keepis a type of lambda function for storing results of reduction in each row. It is declared as
auto keep =[=] __cuda_callable__ ( IndexType rowIdx, const RealType& value ) { ... };
Template Parameters
FetchValueis type returned by the Fetch lambda function.
Parameters
begindefines beginning of the range [begin,end) of rows to be processed.
enddefines ending of the range [begin,end) of rows to be processed.
fetchis an instance of lambda function for data fetch.
reduceis an instance of lambda function for reduction.
keepin an instance of lambda function for storing results.
identityis the identity element for the reduction operation, i.e. element which does not change the result of the reduction.
Example
#include <iostream>
#include <iomanip>
#include <functional>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void reduceRows()
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 0 / 1 . . . . \ -> { 0, 0, 1 }
* 0 | 2 1 . . . | -> { 0, 2, 1 }
* | 3 2 1 . . | -> { 3, 2, 1 }
* | . 3 2 1 . | -> { 3, 2, 1 }
* \ . . 3 2 1 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -2, -1, 0 }.
*/
5, // number of matrix columns
{ -2, -1, 0 }, // diagonals offsets
{ { 0, 0, 1 }, // matrix elements
{ 0, 2, 1 },
{ 3, 2, 1 },
{ 3, 2, 1 },
{ 3, 2, 1 } } );
/***
* Find largest element in each row.
*/
TNL::Containers::Vector< double, Device > rowMax( matrix.getRows() );
/***
* Prepare vector view for lambdas.
*/
auto rowMaxView = rowMax.getView();
/***
* Fetch lambda just returns absolute value of matrix elements.
*/
auto fetch = [=] __cuda_callable__ ( int rowIdx, int columnIdx, const double& value ) -> double {
return TNL::abs( value );
};
/***
* Reduce lambda return maximum of given values.
*/
auto reduce = [=] __cuda_callable__ ( const double& a, const double& b ) -> double {
return TNL::max( a, b );
};
/***
* Keep lambda store the largest value in each row to the vector rowMax.
*/
auto keep = [=] __cuda_callable__ ( int rowIdx, const double& value ) mutable {
rowMaxView[ rowIdx ] = value;
};
/***
* Compute the largest values in each row.
*/
matrix.reduceRows( 0, matrix.getRows(), fetch, reduce, keep, std::numeric_limits< double >::lowest() );
std::cout << "The matrix reads as: " << std::endl << matrix << std::endl;
std::cout << "Max. elements in rows are: " << rowMax << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Rows reduction on host:" << std::endl;
reduceRows< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Rows reduction on CUDA device:" << std::endl;
reduceRows< TNL::Devices::Cuda >();
#endif
}
Output
Rows reduction on host:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1
Max. elements in rows are: [ 1, 2, 3, 3, 3 ]
Rows reduction on CUDA device:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 1:3 2:2 3:1
Row: 4 -> 2:3 3:2 4:1
Max. elements in rows are: [ 1, 2, 3, 3, 3 ]

◆ 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::Matrices::Matrix< Real, Device, Index, RealAllocator >.

◆ sequentialForAllRows() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::sequentialForAllRows ( Function &  function)

This method calls sequentialForRows for all matrix rows.

See MultidiagonalMatrix::sequentialForAllRows.

Template Parameters
Functionis a type of lambda function that will operate on matrix elements.
Parameters
functionis an instance of the lambda function to be called in each row.

◆ sequentialForAllRows() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::sequentialForAllRows ( Function &  function) const

This method calls sequentialForRows for all matrix rows (for constant instances).

See MultidiagonalMatrix::sequentialForRows.

Template Parameters
Functionis a type of lambda function that will operate on matrix elements.
Parameters
functionis an instance of the lambda function to be called in each row.

◆ sequentialForRows() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::sequentialForRows ( IndexType  begin,
IndexType  end,
Function &  function 
)

Method for sequential iteration over all matrix rows for non-constant instances.

Template Parameters
Functionis type of lambda function that will operate on matrix elements. It is should have form like
auto function = [=] __cuda_callable__ ( IndexType rowIdx, IndexType columnIdx, IndexType columnIdx_, RealType& value ) {
... };

The column index repeats twice only for compatibility with sparse matrices.

Parameters
begindefines beginning of the range [begin,end) of rows to be processed.
enddefines ending of the range [begin,end) of rows to be processed.
functionis an instance of the lambda function to be called in each row.

◆ sequentialForRows() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::sequentialForRows ( IndexType  begin,
IndexType  end,
Function &  function 
) const

Method for sequential iteration over all matrix rows for constant instances.

Template Parameters
Functionis type of lambda function that will operate on matrix elements. It is should have form like
auto function = [=] ( IndexType rowIdx, IndexType columnIdx, IndexType columnIdx_, const RealType& value ) { ... };

The column index repeats twice only for compatibility with sparse matrices.

Parameters
begindefines beginning of the range [begin,end) of rows to be processed.
enddefines ending of the range [begin,end) of rows to be processed.
functionis an instance of the lambda function to be called in each row.

◆ setDiagonalsOffsets() [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 >::setDiagonalsOffsets ( const Vector &  diagonalsOffsets)

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
diagonalsOffsetsis a vector-like container holding the subdiagonals offsets.

◆ setDiagonalsOffsets() [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 >::setDiagonalsOffsets ( std::initializer_list< ListIndex >  diagonalsOffsets)

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
diagonalsOffsetsis 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 ( IndexType  rows,
IndexType  columns,
const Vector &  diagonalsOffsets 
)

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.
diagonalsOffsetsis vector with diagonals offsets.

◆ setElement()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
__cuda_callable__ void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::setElement ( IndexType  row,
IndexType  column,
const RealType value 
)

Sets element at given row and column to given value.

This method can be called from the host system (CPU) no matter where the matrix is allocated. If the matrix is allocated on GPU this method can be called even from device kernels. If the matrix is allocated in GPU device this method is called from CPU, it transfers values of each matrix element separately and so the performance is very low. For higher performance see. MultidiagonalMatrix::getRow or MultidiagonalMatrix::forElements and MultidiagonalMatrix::forAllElements. The call may fail if the matrix row capacity is exhausted.

Parameters
rowis row index of the element.
columnis columns index of the element.
valueis the value the element will be set to.
Example
#include <iostream>
#include <TNL/Algorithms/parallelFor.h>
#include <TNL/Matrices/MultidiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
#include <TNL/Pointers/SharedPointer.h>
#include <TNL/Pointers/SmartPointersRegister.h>
template< typename Device >
{
const int matrixSize( 5 );
auto diagonalsOffsets = { -1, 0, 1 }; // offsets of tridiagonal matrix
TNL::Pointers::SharedPointer< Matrix > matrix( matrixSize, matrixSize, diagonalsOffsets );
for( int i = 0; i < 5; i++ )
matrix->setElement( i, i, i );
std::cout << "Matrix set from the host:" << std::endl;
std::cout << *matrix << std::endl;
Matrix* matrix_device = &matrix.template modifyData< Device >();
auto f = [=] __cuda_callable__ ( int i ) mutable {
if( i > 0 )
matrix_device->setElement( i, i - 1, 1.0 );
matrix_device->setElement( i, i, -i );
if( i < matrixSize - 1 )
matrix_device->setElement( i, i + 1, 1.0 );
};
/***
* For the case when Device is CUDA device we need to synchronize smart
* pointers. To avoid this you may use MultidiagonalMatrixView. See
* MultidiagonalMatrixView::getRow example for details.
*/
TNL::Pointers::synchronizeSmartPointersOnDevice< Device >();
TNL::Algorithms::parallelFor< Device >( 0, matrixSize, f );
std::cout << "Matrix set from its native device:" << std::endl;
std::cout << *matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Set elements on host:" << std::endl;
setElements< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Set elements on CUDA device:" << std::endl;
setElements< TNL::Devices::Cuda >();
#endif
}
Base class for other matrix types.
Definition: Matrix.h:38
Output
Set elements on host:
Matrix set from the host:
Row: 0 ->
Row: 1 -> 1:1
Row: 2 -> 2:2
Row: 3 -> 3:3
Row: 4 -> 4:4
Matrix set from its native device:
Row: 0 -> 1:1
Row: 1 -> 0:1 1:-1 2:1
Row: 2 -> 1:1 2:-2 3:1
Row: 3 -> 2:1 3:-3 4:1
Row: 4 -> 3:1 4:-4
Set elements on CUDA device:
Matrix set from the host:
Row: 0 ->
Row: 1 -> 1:1
Row: 2 -> 2:2
Row: 3 -> 3:3
Row: 4 -> 4:4
Matrix set from its native device:
Row: 0 -> 1:1
Row: 1 -> 0:1 1:-1 2:1
Row: 2 -> 1:1 2:-2 3:1
Row: 3 -> 2:1 3:-3 4:1
Row: 4 -> 3:1 4:-4

◆ 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 }.
*/
matrixSize, // number of matrix rows
matrixSize, // number of matrix columns
{ -3, -1, 0, 1, 3 } ); // matrix diagonals offsets
matrix.setElements( {
/***
* 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 }
*
*/
{ 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 }
} );
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

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

◆ setValue()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::setValue ( const RealType value)

Set all matrix elements to given value.

Parameters
valueis the new value of all matrix elements.

◆ vectorProduct()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator , typename IndexAllocator >
template<typename InVector , typename OutVector >
void TNL::Matrices::MultidiagonalMatrix< Real, Device, Index, Organization, RealAllocator, IndexAllocator >::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.

More precisely, it computes:

outVector = matrixMultiplicator * ( * this ) * inVector + outVectorMultiplicator * outVector
Template Parameters
InVectoris type of input vector. It can be TNL::Containers::Vector, TNL::Containers::VectorView, TNL::Containers::Array, TNL::Containers::ArrayView, or similar container.
OutVectoris type of output vector. It can be TNL::Containers::Vector, TNL::Containers::VectorView, TNL::Containers::Array, TNL::Containers::ArrayView, or similar container.
Parameters
inVectoris input vector.
outVectoris output vector.
matrixMultiplicatoris a factor by which the matrix is multiplied. It is one by default.
outVectorMultiplicatoris a factor by which the outVector is multiplied before added to the result of matrix-vector product. It is zero by default.
beginis the beginning of the rows range for which the vector product is computed. It is zero by default.
endis the end of the rows range for which the vector product is computed. It is number if the matrix rows by default.

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