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::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator > Class Template Reference

Implementation of sparse matrix, i.e. matrix storing only non-zero elements. More...

#include <TNL/Matrices/SparseMatrix.h>

Inheritance diagram for TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >:
Inheritance graph
[legend]
Collaboration diagram for TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >:
Collaboration graph
[legend]

Public Types

using BaseType = Matrix< Real, Device, Index, RealAllocator >
 
using ColumnsIndexesVectorType = Containers::Vector< typename TNL::copy_const< Index >::template from< Real >::type, Device, Index, IndexAllocator >
 
using ColumnsIndexesViewType = typename ColumnsIndexesVectorType::ViewType
 
using ComputeRealType = ComputeReal
 
using ConstColumnsIndexesViewType = typename ColumnsIndexesViewType::ConstViewType
 
using ConstRowsCapacitiesView = typename RowsCapacitiesView::ConstViewType
 
using ConstRowView = typename ViewType::ConstRowView
 Type for accessing constant matrix rows.
 
using ConstValuesViewType = typename ValuesViewType::ConstViewType
 
using ConstViewType = SparseMatrixView< std::add_const_t< Real >, Device, Index, MatrixType, SegmentsViewTemplate, ComputeRealType >
 Matrix view type for constant instances. More...
 
using DeviceType = Device
 The device where the matrix is allocated.
 
using IndexAllocatorType = IndexAllocator
 The allocator for matrix elements column indexes.
 
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 RowsCapacitiesType = Containers::Vector< std::remove_const_t< Index >, Device, Index, IndexAllocator >
 
using RowsCapacitiesView = Containers::VectorView< std::remove_const_t< Index >, Device, Index >
 
using RowView = typename ViewType::RowView
 Type for accessing matrix rows.
 
template<typename Device_ , typename Index_ , typename IndexAllocator_ >
using SegmentsTemplate = Segments< Device_, Index_, IndexAllocator_ >
 Templated type of segments, i.e. sparse matrix format.
 
using SegmentsType = Segments< Device, Index, IndexAllocator >
 Type of segments used by this matrix. It represents the sparse matrix format.
 
template<typename Device_ , typename Index_ >
using SegmentsViewTemplate = typename SegmentsType::template ViewTemplate< Device_, Index >
 Templated view type of segments, i.e. sparse matrix format.
 
using SegmentsViewType = typename SegmentsType::ViewType
 Type of segments view used by the related matrix view. It represents the sparse matrix format.
 
template<typename _Real = Real, typename _Device = Device, typename _Index = Index, typename _MatrixType = MatrixType, template< typename, typename, typename > class _Segments = Segments, typename _ComputeReal = ComputeReal, typename _RealAllocator = typename Allocators::Default< _Device >::template Allocator< _Real >, typename _IndexAllocator = typename Allocators::Default< _Device >::template Allocator< _Index >>
using Self = SparseMatrix< _Real, _Device, _Index, _MatrixType, _Segments, _ComputeReal, _RealAllocator, _IndexAllocator >
 Helper type for getting self type or its modifications.
 
using ValuesVectorType = typename Matrix< Real, Device, Index, RealAllocator >::ValuesType
 
using ValuesViewType = typename ValuesVectorType::ViewType
 
using ViewType = SparseMatrixView< Real, Device, Index, MatrixType, SegmentsViewTemplate, ComputeRealType >
 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

 SparseMatrix (const RealAllocatorType &realAllocator=RealAllocatorType(), const IndexAllocatorType &indexAllocator=IndexAllocatorType())
 Constructor only with values and column indexes allocators. More...
 
template<typename RowCapacitiesVector , std::enable_if_t< TNL::IsArrayType< RowCapacitiesVector >::value, int > = 0>
 SparseMatrix (const RowCapacitiesVector &rowCapacities, IndexType columns, const RealAllocatorType &realAllocator=RealAllocatorType(), const IndexAllocatorType &indexAllocator=IndexAllocatorType())
 Constructor with matrix rows capacities given as a vector and number of columns. More...
 
 SparseMatrix (const SparseMatrix &matrix)
 Copy constructor. More...
 
template<typename ListIndex >
 SparseMatrix (const std::initializer_list< ListIndex > &rowCapacities, IndexType columns, const RealAllocatorType &realAllocator=RealAllocatorType(), const IndexAllocatorType &indexAllocator=IndexAllocatorType())
 Constructor with matrix rows capacities and number of columns. More...
 
template<typename Index_t , std::enable_if_t< std::is_integral< Index_t >::value, int > = 0>
 SparseMatrix (Index_t rows, Index_t columns, const RealAllocatorType &realAllocator=RealAllocatorType(), const IndexAllocatorType &indexAllocator=IndexAllocatorType())
 Constructor with matrix dimensions. More...
 
 SparseMatrix (IndexType rows, IndexType columns, const std::initializer_list< std::tuple< IndexType, IndexType, RealType > > &data, const RealAllocatorType &realAllocator=RealAllocatorType(), const IndexAllocatorType &indexAllocator=IndexAllocatorType())
 Constructor with matrix dimensions and data in initializer list. More...
 
template<typename MapIndex , typename MapValue >
 SparseMatrix (IndexType rows, IndexType columns, const std::map< std::pair< MapIndex, MapIndex >, MapValue > &map, const RealAllocatorType &realAllocator=RealAllocatorType(), const IndexAllocatorType &indexAllocator=IndexAllocatorType())
 Constructor with matrix dimensions and data in std::map. More...
 
 SparseMatrix (SparseMatrix &&matrix) noexcept=default
 Move constructor. More...
 
__cuda_callable__ void addElement (IndexType row, IndexType column, const RealType &value, const RealType &thisElementMultiplicator)
 Add element at given row and column to given value. More...
 
template<typename Function >
void forAllElements (Function &&function)
 Method for parallel iteration over all matrix elements for non-constant instances. More...
 
template<typename Function >
void forAllElements (Function &&function) const
 Method for parallel iteration over all matrix elements 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 parallel iteration over all matrix elements of given rows for non-constant instances. More...
 
template<typename Function >
void forElements (IndexType begin, IndexType end, Function &&function) const
 Method for parallel iteration over matrix elements of given 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...
 
ColumnsIndexesVectorTypegetColumnIndexes ()
 Getter of column indexes for nonconstant instances. More...
 
const ColumnsIndexesVectorTypegetColumnIndexes () const
 Getter of column indexes 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 sparse matrix. 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...
 
IndexType getNonzeroElementsCount () const override
 Returns number of non-zero matrix elements. More...
 
__cuda_callable__ IndexType getPaddingIndex () const
 Returns a padding index value. 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...
 
__cuda_callable__ IndexType getRowCapacity (IndexType row) const
 Returns capacity of given matrix row. More...
 
SegmentsTypegetSegments ()
 Getter of segments for non-constant instances. More...
 
const SegmentsTypegetSegments () const
 Getter of segments for constant instances. More...
 
std::string getSerializationTypeVirtual () const override
 Returns string with serialization type. More...
 
template<typename Real2 , typename Index2 , template< typename, typename, typename > class Segments2>
void getTransposition (const SparseMatrix< Real2, Device, Index2, MatrixType, Segments2 > &matrix, const ComputeRealType &matrixMultiplicator=1.0)
 
ViewType getView ()
 Returns a modifiable view of the sparse 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 Matrix >
bool operator!= (const Matrix &matrix) const
 Comparison operator with another arbitrary matrix type. More...
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization, typename RealAllocator_ >
SparseMatrixoperator= (const DenseMatrix< Real_, Device_, Index_, Organization, RealAllocator_ > &matrix)
 Assignment of dense matrix. More...
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization, typename RealAllocator_ >
SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator > & operator= (const DenseMatrix< Real_, Device_, Index_, Organization, RealAllocator_ > &matrix)
 
template<typename RHSMatrix >
SparseMatrixoperator= (const RHSMatrix &matrix)
 Assignment of any matrix type other then this and dense. More...
 
template<typename RHSMatrix >
SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator > & operator= (const RHSMatrix &matrix)
 
SparseMatrixoperator= (const SparseMatrix &matrix)
 Copy-assignment of exactly the same matrix type. More...
 
SparseMatrixoperator= (SparseMatrix &&matrix) noexcept(false)
 Move-assignment of exactly the same matrix type. More...
 
template<typename Matrix >
bool operator== (const Matrix &matrix) const
 Comparison operator with another arbitrary matrix type. 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, const 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, const 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, const 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, const Reduce &reduce, Keep &keep, const FetchReal &identity) const
 Method for performing general reduction on matrix rows for constant instances. More...
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchValue >
void reduceRows (IndexType begin, IndexType end, Fetch &fetch, const Reduce &reduce, Keep &keep, const FetchValue &identity)
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchValue >
void reduceRows (IndexType begin, IndexType end, Fetch &fetch, const Reduce &reduce, Keep &keep, const FetchValue &identity) const
 
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...
 
virtual void setColumnsWithoutReset (IndexType columns)
 Set number of columns of this matrix. More...
 
void setDimensions (IndexType rows, IndexType columns) override
 Set number of rows and columns of this matrix. More...
 
__cuda_callable__ void setElement (IndexType row, IndexType column, const RealType &value)
 Sets element at given row and column to given value. More...
 
void setElements (const std::initializer_list< std::tuple< IndexType, IndexType, RealType > > &data)
 This method sets the sparse matrix elements from initializer list. More...
 
template<typename MapIndex , typename MapValue >
void setElements (const std::map< std::pair< MapIndex, MapIndex >, MapValue > &map)
 This method sets the sparse matrix elements from std::map. More...
 
template<typename Matrix >
void setLike (const Matrix &matrix)
 Set the number of matrix rows and columns by the given matrix. More...
 
template<typename Matrix_ >
void setLike (const Matrix_ &matrix)
 
template<typename RowsCapacitiesVector >
void setRowCapacities (const RowsCapacitiesVector &rowCapacities)
 Allocates memory for non-zero matrix elements. More...
 
template<typename InVector , typename OutVector >
void vectorProduct (const InVector &inVector, OutVector &outVector, ComputeRealType matrixMultiplicator=1.0, ComputeRealType 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 std::string getSerializationType ()
 Returns string with serialization type. More...
 
static constexpr bool isBinary ()
 Test of binary matrix type. More...
 
static constexpr bool isSymmetric ()
 Test of symmetric matrix type. More...
 
- Static Public Member Functions inherited from TNL::Object
static std::string getSerializationType ()
 Static serialization type getter. More...
 

Protected Attributes

ColumnsIndexesVectorType columnIndexes
 
IndexAllocator indexAllocator
 
SegmentsType segments
 
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, typename MatrixType = GeneralMatrix, template< typename Device_, typename Index_, typename IndexAllocator_ > class Segments = Algorithms::Segments::CSRDefault, typename ComputeReal = typename ChooseSparseMatrixComputeReal< Real, Index >::type, typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
class TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >

Implementation of sparse matrix, i.e. matrix storing only non-zero elements.

Template Parameters
Realis a type of matrix elements. If Real equals bool the matrix is treated as binary and so the matrix elements values are not stored in the memory since we need to remember only coordinates of non-zero elements( which equal one).
Deviceis a device where the matrix is allocated.
Indexis a type for indexing of the matrix elements.
MatrixTypespecifies a symmetry of matrix. See MatrixType. Symmetric matrices store only lower part of the matrix and its diagonal. The upper part is reconstructed on the fly. GeneralMatrix with no symmetry is used by default.
Segmentsis a structure representing the sparse matrix format. Depending on the pattern of the non-zero elements different matrix formats can perform differently especially on GPUs. By default Algorithms::Segments::CSR format is used. See also Algorithms::Segments::Ellpack, Algorithms::Segments::SlicedEllpack, Algorithms::Segments::ChunkedEllpack, and Algorithms::Segments::BiEllpack.
ComputeRealis the same as Real mostly but for binary matrices it is set to Index type. This can be changed by the user, of course.
RealAllocatoris allocator for the matrix elements values.
IndexAllocatoris allocator for the matrix elements column indexes.

Member Typedef Documentation

◆ ConstViewType

template<typename Real = double, typename Device = Devices::Host, typename Index = int, typename MatrixType = GeneralMatrix, template< typename Device_, typename Index_, typename IndexAllocator_ > class Segments = Algorithms::Segments::CSRDefault, typename ComputeReal = typename ChooseSparseMatrixComputeReal< Real, Index >::type, typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
using TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::ConstViewType = SparseMatrixView< std::add_const_t< Real >, Device, Index, MatrixType, SegmentsViewTemplate, ComputeRealType >

Matrix view type for constant instances.

See SparseMatrixView.

◆ ViewType

template<typename Real = double, typename Device = Devices::Host, typename Index = int, typename MatrixType = GeneralMatrix, template< typename Device_, typename Index_, typename IndexAllocator_ > class Segments = Algorithms::Segments::CSRDefault, typename ComputeReal = typename ChooseSparseMatrixComputeReal< Real, Index >::type, typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
using TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::ViewType = SparseMatrixView< Real, Device, Index, MatrixType, SegmentsViewTemplate, ComputeRealType >

Type of related matrix view.

See SparseMatrixView.

Constructor & Destructor Documentation

◆ SparseMatrix() [1/8]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::SparseMatrix ( const RealAllocatorType realAllocator = RealAllocatorType(),
const IndexAllocatorType indexAllocator = IndexAllocatorType() 
)

Constructor only with values and column indexes allocators.

Parameters
realAllocatoris used for allocation of matrix elements values.
indexAllocatoris used for allocation of matrix elements column indexes.

◆ SparseMatrix() [2/8]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::SparseMatrix ( const SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator > &  matrix)
explicit

Copy constructor.

Parameters
matrixis the source matrix

◆ SparseMatrix() [3/8]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, typename MatrixType = GeneralMatrix, template< typename Device_, typename Index_, typename IndexAllocator_ > class Segments = Algorithms::Segments::CSRDefault, typename ComputeReal = typename ChooseSparseMatrixComputeReal< Real, Index >::type, typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::SparseMatrix ( SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator > &&  matrix)
defaultnoexcept

Move constructor.

Parameters
matrixis the source matrix

◆ SparseMatrix() [4/8]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Index_t , std::enable_if_t< std::is_integral< Index_t >::value, int > >
TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::SparseMatrix ( Index_t  rows,
Index_t  columns,
const RealAllocatorType realAllocator = RealAllocatorType(),
const IndexAllocatorType indexAllocator = IndexAllocatorType() 
)

Constructor with matrix dimensions.

Parameters
rowsis number of matrix rows.
columnsis number of matrix columns.
realAllocatoris used for allocation of matrix elements values.
indexAllocatoris used for allocation of matrix elements column indexes.

◆ SparseMatrix() [5/8]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename ListIndex >
TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::SparseMatrix ( const std::initializer_list< ListIndex > &  rowCapacities,
IndexType  columns,
const RealAllocatorType realAllocator = RealAllocatorType(),
const IndexAllocatorType indexAllocator = IndexAllocatorType() 
)
explicit

Constructor with matrix rows capacities and number of columns.

The number of matrix rows is given by the size of rowCapacities list.

Template Parameters
ListIndexis the initializer list values type.
Parameters
rowCapacitiesis a list telling how many matrix elements must be allocated in each row.
columnsis the number of matrix columns.
realAllocatoris used for allocation of matrix elements values.
indexAllocatoris used for allocation of matrix elements column indexes.
Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void initializerListExample()
{
{ 1, 2, 3, 4, 5 }, // row capacities
6 }; // number of matrix columns
for( int row = 0; row < matrix.getRows(); row++ )
for( int column = 0; column <= row; column++ )
matrix.setElement( row, column, row - column + 1 );
std::cout << "General sparse matrix: " << std::endl << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrices on CPU ... " << std::endl;
initializerListExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrices on CUDA GPU ... " << std::endl;
initializerListExample< TNL::Devices::Cuda >();
#endif
}
Implementation of sparse matrix, i.e. matrix storing only non-zero elements.
Definition: SparseMatrix.h:57
__cuda_callable__ void setElement(IndexType row, IndexType column, const RealType &value)
Sets element at given row and column to given value.
Definition: SparseMatrix.hpp:502
T endl(T... args)
Output
Creating matrices on CPU ...
General sparse matrix:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 0:4 1:3 2:2 3:1
Row: 4 -> 0:5 1:4 2:3 3:2 4:1
Creating matrices on CUDA GPU ...
General sparse matrix:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 0:4 1:3 2:2 3:1
Row: 4 -> 0:5 1:4 2:3 3:2 4:1

◆ SparseMatrix() [6/8]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename RowCapacitiesVector , std::enable_if_t< TNL::IsArrayType< RowCapacitiesVector >::value, int > >
TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::SparseMatrix ( const RowCapacitiesVector &  rowCapacities,
IndexType  columns,
const RealAllocatorType realAllocator = RealAllocatorType(),
const IndexAllocatorType indexAllocator = IndexAllocatorType() 
)
explicit

Constructor with matrix rows capacities given as a vector and number of columns.

The number of matrix rows is given by the size of rowCapacities vector.

Template Parameters
RowCapacitiesVectoris the row capacities vector type. Usually it is some of TNL::Containers::Array, TNL::Containers::ArrayView, TNL::Containers::Vector or TNL::Containers::VectorView.
Parameters
rowCapacitiesis a vector telling how many matrix elements must be allocated in each row.
columnsis the number of matrix columns.
realAllocatoris used for allocation of matrix elements values.
indexAllocatoris used for allocation of matrix elements column indexes.
Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Containers/Vector.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void initializerListExample()
{
TNL::Containers::Vector< int, Device > rowCapacities{ 1, 2, 3, 4, 5 };
rowCapacities, // row capacities
6 }; // number of matrix columns
for( int row = 0; row < matrix.getRows(); row++ )
for( int column = 0; column <= row; column++ )
matrix.setElement( row, column, row - column + 1 );
std::cout << "General sparse matrix: " << std::endl << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrices on CPU ... " << std::endl;
initializerListExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrices on CUDA GPU ... " << std::endl;
initializerListExample< TNL::Devices::Cuda >();
#endif
}
Vector extends Array with algebraic operations.
Definition: Vector.h:40
Output
Creating matrices on CPU ...
General sparse matrix:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 0:4 1:3 2:2 3:1
Row: 4 -> 0:5 1:4 2:3 3:2 4:1
Creating matrices on CUDA GPU ...
General sparse matrix:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 0:4 1:3 2:2 3:1
Row: 4 -> 0:5 1:4 2:3 3:2 4:1

◆ SparseMatrix() [7/8]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::SparseMatrix ( IndexType  rows,
IndexType  columns,
const std::initializer_list< std::tuple< IndexType, IndexType, RealType > > &  data,
const RealAllocatorType realAllocator = RealAllocatorType(),
const IndexAllocatorType indexAllocator = IndexAllocatorType() 
)
explicit

Constructor with matrix dimensions and data in initializer list.

The matrix elements values are given as a list data of triples: { { row1, column1, value1 }, { row2, column2, value2 }, ... }.

Parameters
rowsis number of matrix rows.
columnsis number of matrix columns.
datais a list of matrix elements values.
realAllocatoris used for allocation of matrix elements values.
indexAllocatoris used for allocation of matrix elements column indexes.
Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void initializerListExample()
{
5, // number of matrix rows
5, // number of matrix columns
{ // matrix elements definition
{ 0, 0, 2.0 },
{ 1, 0, -1.0 }, { 1, 1, 2.0 }, { 1, 2, -1.0 },
{ 2, 1, -1.0 }, { 2, 2, 2.0 }, { 2, 3, -1.0 },
{ 3, 2, -1.0 }, { 3, 3, 2.0 }, { 3, 4, -1.0 },
{ 4, 4, 2.0 } } );
std::cout << "General sparse matrix: " << std::endl << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrices on CPU ... " << std::endl;
initializerListExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrices on CUDA GPU ... " << std::endl;
initializerListExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrices on CPU ...
General sparse matrix:
Row: 0 -> 0:2
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 -> 4:2
Creating matrices on CUDA GPU ...
General sparse matrix:
Row: 0 -> 0:2
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 -> 4:2

◆ SparseMatrix() [8/8]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename MapIndex , typename MapValue >
TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::SparseMatrix ( IndexType  rows,
IndexType  columns,
const std::map< std::pair< MapIndex, MapIndex >, MapValue > &  map,
const RealAllocatorType realAllocator = RealAllocatorType(),
const IndexAllocatorType indexAllocator = IndexAllocatorType() 
)
explicit

Constructor with matrix dimensions and data in std::map.

The matrix elements values are given as a map data where keys are std::pair of matrix coordinates ( {row, column} ) and value is the matrix element value.

Template Parameters
MapIndexis a type for indexing rows and columns.
MapValueis a type for matrix elements values in the map.
Parameters
rowsis number of matrix rows.
columnsis number of matrix columns.
mapis std::map containing matrix elements.
realAllocatoris used for allocation of matrix elements values.
indexAllocatoris used for allocation of matrix elements column indexes.
Example
#include <iostream>
#include <map>
#include <utility>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void initializerListExample()
{
map.insert( std::make_pair( std::make_pair( 0, 0 ), 2.0 ) );
map.insert( std::make_pair( std::make_pair( 1, 0 ), -1.0 ) );
map.insert( std::make_pair( std::make_pair( 1, 1 ), 2.0 ) );
map.insert( std::make_pair( std::make_pair( 1, 2 ), -1.0 ) );
map.insert( std::make_pair( std::make_pair( 2, 1 ), -1.0 ) );
map.insert( std::make_pair( std::make_pair( 2, 2 ), 2.0 ) );
map.insert( std::make_pair( std::make_pair( 2, 3 ), -1.0 ) );
map.insert( std::make_pair( std::make_pair( 3, 2 ), -1.0 ) );
map.insert( std::make_pair( std::make_pair( 3, 3 ), 2.0 ) );
map.insert( std::make_pair( std::make_pair( 3, 4 ), -1.0 ) );
map.insert( std::make_pair( std::make_pair( 4, 4 ), 2.0 ) );
std::cout << "General sparse matrix: " << std::endl << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrices on CPU ... " << std::endl;
initializerListExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrices on CUDA GPU ... " << std::endl;
initializerListExample< TNL::Devices::Cuda >();
#endif
}
T insert(T... args)
T make_pair(T... args)
Output
Creating matrices on CPU ...
General sparse matrix:
Row: 0 -> 0:2
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 -> 4:2
Creating matrices on CUDA GPU ...
General sparse matrix:
Row: 0 -> 0:2
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 -> 4:2

Member Function Documentation

◆ addElement()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
__cuda_callable__ void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::addElement ( IndexType  row,
IndexType  column,
const RealType value,
const RealType thisElementMultiplicator 
)

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. SparseMatrix::getRow or SparseMatrix::forElements and SparseMatrix::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/SparseMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void addElements()
{
TNL::Matrices::SparseMatrix< double, Device > matrix( { 5, 5, 5, 5, 5 }, 5 );
for( int i = 0; i < 5; i++ )
matrix.setElement( i, i, i );
std::cout << "Initial matrix is: " << std::endl << matrix << std::endl;
for( int i = 0; i < 5; i++ )
for( int j = 0; j < 5; j++ )
matrix.addElement( i, j, 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
}
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 2:1 3:1 4:1
Row: 1 -> 0:1 1:6 2:1 3:1 4:1
Row: 2 -> 0:1 1:1 2:11 3:1 4:1
Row: 3 -> 0:1 1:1 2:1 3:16 4:1
Row: 4 -> 0:1 1:1 2:1 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 2:1 3:1 4:1
Row: 1 -> 0:1 1:6 2:1 3:1 4:1
Row: 2 -> 0:1 1:1 2:11 3:1 4:1
Row: 3 -> 0:1 1:1 2:1 3:16 4:1
Row: 4 -> 0:1 1:1 2:1 3:1 4:21

◆ forAllElements() [1/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::forAllElements ( Function &&  function)

Method for parallel iteration over all matrix elements for non-constant instances.

See SparseMatrix::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 for each matrix element.
Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forElementsExample()
{
TNL::Matrices::SparseMatrix< double, Device > matrix( { 1, 2, 3, 4, 5 }, 5 );
auto f = [=] __cuda_callable__ ( int rowIdx, int localIdx, int& columnIdx, double& value ) {
if( rowIdx >= localIdx ) // This is important, some matrix formats may allocate more matrix elements
// than we requested. These padding elements are processed here as well.
{
columnIdx = localIdx;
value = rowIdx + localIdx + 1;
}
};
matrix.forElements( 0, matrix.getRows(), f );
std::cout << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forElementsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forElementsExample< TNL::Devices::Cuda >();
#endif
}
#define __cuda_callable__
Definition: CudaCallable.h:22
Output
Creating matrix on host:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:3
Row: 2 -> 0:3 1:4 2:5
Row: 3 -> 0:4 1:5 2:6 3:7
Row: 4 -> 0:5 1:6 2:7 3:8 4:9
Creating matrix on CUDA device:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:3
Row: 2 -> 0:3 1:4 2:5
Row: 3 -> 0:4 1:5 2:6 3:7
Row: 4 -> 0:5 1:6 2:7 3:8 4:9

◆ forAllElements() [2/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::forAllElements ( Function &&  function) const

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

See SparseMatrix::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 for each matrix element.
Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forElementsExample()
{
TNL::Matrices::SparseMatrix< double, Device > matrix( { 1, 2, 3, 4, 5 }, 5 );
auto f = [=] __cuda_callable__ ( int rowIdx, int localIdx, int& columnIdx, double& value ) {
if( rowIdx >= localIdx ) // This is important, some matrix formats may allocate more matrix elements
// than we requested. These padding elements are processed here as well.
{
columnIdx = localIdx;
value = rowIdx + localIdx + 1;
}
};
matrix.forElements( 0, matrix.getRows(), f );
std::cout << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forElementsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forElementsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:3
Row: 2 -> 0:3 1:4 2:5
Row: 3 -> 0:4 1:5 2:6 3:7
Row: 4 -> 0:5 1:6 2:7 3:8 4:9
Creating matrix on CUDA device:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:3
Row: 2 -> 0:3 1:4 2:5
Row: 3 -> 0:4 1:5 2:6 3:7
Row: 4 -> 0:5 1:6 2:7 3:8 4:9

◆ forAllRows() [1/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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 SparseMatrix::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: SparseMatrix.h:167

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

Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forRowsExample()
{
/***
* Set the following matrix (dots represent zero matrix elements):
*
* / 2 . . . . \
* | 1 2 1 . . |
* | . 1 2 1. . |
* | . . 1 2 1 |
* \ . . . . 2 /
*/
const int size( 5 );
MatrixType matrix( { 1, 3, 3, 3, 1 }, size );
using RowView = typename MatrixType::RowView;
/***
* Set the matrix elements.
*/
auto f = [=] __cuda_callable__ ( RowView& row ) mutable {
const int rowIdx = row.getRowIndex();
if( rowIdx == 0 )
row.setElement( 0, rowIdx, 2.0 ); // diagonal element
else if( rowIdx == size - 1 )
row.setElement( 0, rowIdx, 2.0 ); // diagonal element
else
{
row.setElement( 0, rowIdx - 1, 1.0 ); // elements below the diagonal
row.setElement( 1, rowIdx, 2.0 ); // diagonal element
row.setElement( 2, rowIdx + 1, 1.0 ); // elements above the diagonal
}
};
matrix.forAllRows( f );
std::cout << matrix << std::endl;
/***
* Divide each matrix row by a sum of all elements in the row - with use of iterators.
*/
matrix.forAllRows( [=] __cuda_callable__ ( RowView& row ) mutable {
double sum( 0.0 );
for( auto element : row )
sum += element.value();
for( auto element: row )
element.value() /= sum;
} );
std::cout << "Divide each matrix row by a sum of all elements in the row ... " << std::endl;
std::cout << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Getting matrix rows on host: " << std::endl;
forRowsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Getting matrix rows on CUDA device: " << std::endl;
forRowsExample< TNL::Devices::Cuda >();
#endif
}
Structure for specifying type of sparse matrix.
Definition: MatrixType.h:21
Output
Getting matrix rows on host:
Row: 0 -> 0:2
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 -> 4:2
Divide each matrix row by a sum of all elements in the row ...
Row: 0 -> 0:1
Row: 1 -> 0:0.25 1:0.5 2:0.25
Row: 2 -> 1:0.25 2:0.5 3:0.25
Row: 3 -> 2:0.25 3:0.5 4:0.25
Row: 4 -> 4:1
Getting matrix rows on CUDA device:
Row: 0 -> 0:2
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 -> 4:2
Divide each matrix row by a sum of all elements in the row ...
Row: 0 -> 0:1
Row: 1 -> 0:0.25 1:0.5 2:0.25
Row: 2 -> 1:0.25 2:0.5 3:0.25
Row: 3 -> 2:0.25 3:0.5 4:0.25
Row: 4 -> 4:1

◆ forAllRows() [2/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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 SparseMatrix::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::SparseMatrix::RowView.

Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forRowsExample()
{
/***
* Set the following matrix (dots represent zero matrix elements):
*
* / 2 . . . . \
* | 1 2 1 . . |
* | . 1 2 1. . |
* | . . 1 2 1 |
* \ . . . . 2 /
*/
const int size( 5 );
MatrixType matrix( { 1, 3, 3, 3, 1 }, size );
using RowView = typename MatrixType::RowView;
/***
* Set the matrix elements.
*/
auto f = [=] __cuda_callable__ ( RowView& row ) mutable {
const int rowIdx = row.getRowIndex();
if( rowIdx == 0 )
row.setElement( 0, rowIdx, 2.0 ); // diagonal element
else if( rowIdx == size - 1 )
row.setElement( 0, rowIdx, 2.0 ); // diagonal element
else
{
row.setElement( 0, rowIdx - 1, 1.0 ); // elements below the diagonal
row.setElement( 1, rowIdx, 2.0 ); // diagonal element
row.setElement( 2, rowIdx + 1, 1.0 ); // elements above the diagonal
}
};
matrix.forAllRows( f );
std::cout << matrix << std::endl;
/***
* Divide each matrix row by a sum of all elements in the row - with use of iterators.
*/
matrix.forAllRows( [=] __cuda_callable__ ( RowView& row ) mutable {
double sum( 0.0 );
for( auto element : row )
sum += element.value();
for( auto element: row )
element.value() /= sum;
} );
std::cout << "Divide each matrix row by a sum of all elements in the row ... " << std::endl;
std::cout << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Getting matrix rows on host: " << std::endl;
forRowsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Getting matrix rows on CUDA device: " << std::endl;
forRowsExample< TNL::Devices::Cuda >();
#endif
}
Output
Getting matrix rows on host:
Row: 0 -> 0:2
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 -> 4:2
Divide each matrix row by a sum of all elements in the row ...
Row: 0 -> 0:1
Row: 1 -> 0:0.25 1:0.5 2:0.25
Row: 2 -> 1:0.25 2:0.5 3:0.25
Row: 3 -> 2:0.25 3:0.5 4:0.25
Row: 4 -> 4:1
Getting matrix rows on CUDA device:
Row: 0 -> 0:2
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 -> 4:2
Divide each matrix row by a sum of all elements in the row ...
Row: 0 -> 0:1
Row: 1 -> 0:0.25 1:0.5 2:0.25
Row: 2 -> 1:0.25 2:0.5 3:0.25
Row: 3 -> 2:0.25 3:0.5 4:0.25
Row: 4 -> 4:1

◆ forElements() [1/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::forElements ( IndexType  begin,
IndexType  end,
Function &&  function 
)

Method for parallel iteration over all matrix elements of given rows for non-constant instances.

Template Parameters
Functionis type of lambda function that will operate on matrix elements.
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 element of given rows.

The lambda function function should be declared like follows:

auto function = [] __cuda_callable__ ( IndexType rowIdx, IndexType localIdx, IndexType& columnIdx, const RealType& value )
mutable { ... }
Index IndexType
The type used for matrix elements indexing.
Definition: SparseMatrix.h:115
Definition: Real.h:17

The localIdx parameter is a rank of the non-zero element in given row.

Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forElementsExample()
{
TNL::Matrices::SparseMatrix< double, Device > matrix( { 1, 2, 3, 4, 5 }, 5 );
auto f = [=] __cuda_callable__ ( int rowIdx, int localIdx, int& columnIdx, double& value ) {
if( rowIdx >= localIdx ) // This is important, some matrix formats may allocate more matrix elements
// than we requested. These padding elements are processed here as well.
{
columnIdx = localIdx;
value = rowIdx + localIdx + 1;
}
};
matrix.forElements( 0, matrix.getRows(), f );
std::cout << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forElementsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forElementsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:3
Row: 2 -> 0:3 1:4 2:5
Row: 3 -> 0:4 1:5 2:6 3:7
Row: 4 -> 0:5 1:6 2:7 3:8 4:9
Creating matrix on CUDA device:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:3
Row: 2 -> 0:3 1:4 2:5
Row: 3 -> 0:4 1:5 2:6 3:7
Row: 4 -> 0:5 1:6 2:7 3:8 4:9

◆ forElements() [2/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::forElements ( IndexType  begin,
IndexType  end,
Function &&  function 
) const

Method for parallel iteration over matrix elements of given rows for constant instances.

Template Parameters
Functionis type of lambda function that will operate on matrix elements.
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 element of given rows.

The lambda function function should be declared like follows:

auto function = [] __cuda_callable__ ( IndexType rowIdx, IndexType localIdx, IndexType columnIdx, const RealType& value )
{ ... };

The localIdx parameter is a rank of the non-zero element in given row.

Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forElementsExample()
{
TNL::Matrices::SparseMatrix< double, Device > matrix( { 1, 2, 3, 4, 5 }, 5 );
auto f = [=] __cuda_callable__ ( int rowIdx, int localIdx, int& columnIdx, double& value ) {
if( rowIdx >= localIdx ) // This is important, some matrix formats may allocate more matrix elements
// than we requested. These padding elements are processed here as well.
{
columnIdx = localIdx;
value = rowIdx + localIdx + 1;
}
};
matrix.forElements( 0, matrix.getRows(), f );
std::cout << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forElementsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forElementsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:3
Row: 2 -> 0:3 1:4 2:5
Row: 3 -> 0:4 1:5 2:6 3:7
Row: 4 -> 0:5 1:6 2:7 3:8 4:9
Creating matrix on CUDA device:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:3
Row: 2 -> 0:3 1:4 2:5
Row: 3 -> 0:4 1:5 2:6 3:7
Row: 4 -> 0:5 1:6 2:7 3:8 4:9

◆ forRows() [1/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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 SparseMatrix::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::SparseMatrix::RowView.

Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forRowsExample()
{
/***
* Set the following matrix (dots represent zero matrix elements):
*
* / 2 . . . . \
* | 1 2 1 . . |
* | . 1 2 1. . |
* | . . 1 2 1 |
* \ . . . . 2 /
*/
const int size( 5 );
MatrixType matrix( { 1, 3, 3, 3, 1 }, size );
using RowView = typename MatrixType::RowView;
/***
* Set the matrix elements.
*/
auto f = [=] __cuda_callable__ ( RowView& row ) mutable {
const int rowIdx = row.getRowIndex();
if( rowIdx == 0 )
row.setElement( 0, rowIdx, 2.0 ); // diagonal element
else if( rowIdx == size - 1 )
row.setElement( 0, rowIdx, 2.0 ); // diagonal element
else
{
row.setElement( 0, rowIdx - 1, 1.0 ); // elements below the diagonal
row.setElement( 1, rowIdx, 2.0 ); // diagonal element
row.setElement( 2, rowIdx + 1, 1.0 ); // elements above the diagonal
}
};
matrix.forAllRows( f );
std::cout << matrix << std::endl;
/***
* Divide each matrix row by a sum of all elements in the row - with use of iterators.
*/
matrix.forAllRows( [=] __cuda_callable__ ( RowView& row ) mutable {
double sum( 0.0 );
for( auto element : row )
sum += element.value();
for( auto element: row )
element.value() /= sum;
} );
std::cout << "Divide each matrix row by a sum of all elements in the row ... " << std::endl;
std::cout << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Getting matrix rows on host: " << std::endl;
forRowsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Getting matrix rows on CUDA device: " << std::endl;
forRowsExample< TNL::Devices::Cuda >();
#endif
}
Output
Getting matrix rows on host:
Row: 0 -> 0:2
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 -> 4:2
Divide each matrix row by a sum of all elements in the row ...
Row: 0 -> 0:1
Row: 1 -> 0:0.25 1:0.5 2:0.25
Row: 2 -> 1:0.25 2:0.5 3:0.25
Row: 3 -> 2:0.25 3:0.5 4:0.25
Row: 4 -> 4:1
Getting matrix rows on CUDA device:
Row: 0 -> 0:2
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 -> 4:2
Divide each matrix row by a sum of all elements in the row ...
Row: 0 -> 0:1
Row: 1 -> 0:0.25 1:0.5 2:0.25
Row: 2 -> 1:0.25 2:0.5 3:0.25
Row: 3 -> 2:0.25 3:0.5 4:0.25
Row: 4 -> 4:1

◆ forRows() [2/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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 SparseMatrix::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::SparseMatrix::RowView.

Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void forRowsExample()
{
/***
* Set the following matrix (dots represent zero matrix elements):
*
* / 2 . . . . \
* | 1 2 1 . . |
* | . 1 2 1. . |
* | . . 1 2 1 |
* \ . . . . 2 /
*/
const int size( 5 );
MatrixType matrix( { 1, 3, 3, 3, 1 }, size );
using RowView = typename MatrixType::RowView;
/***
* Set the matrix elements.
*/
auto f = [=] __cuda_callable__ ( RowView& row ) mutable {
const int rowIdx = row.getRowIndex();
if( rowIdx == 0 )
row.setElement( 0, rowIdx, 2.0 ); // diagonal element
else if( rowIdx == size - 1 )
row.setElement( 0, rowIdx, 2.0 ); // diagonal element
else
{
row.setElement( 0, rowIdx - 1, 1.0 ); // elements below the diagonal
row.setElement( 1, rowIdx, 2.0 ); // diagonal element
row.setElement( 2, rowIdx + 1, 1.0 ); // elements above the diagonal
}
};
matrix.forAllRows( f );
std::cout << matrix << std::endl;
/***
* Divide each matrix row by a sum of all elements in the row - with use of iterators.
*/
matrix.forAllRows( [=] __cuda_callable__ ( RowView& row ) mutable {
double sum( 0.0 );
for( auto element : row )
sum += element.value();
for( auto element: row )
element.value() /= sum;
} );
std::cout << "Divide each matrix row by a sum of all elements in the row ... " << std::endl;
std::cout << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Getting matrix rows on host: " << std::endl;
forRowsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Getting matrix rows on CUDA device: " << std::endl;
forRowsExample< TNL::Devices::Cuda >();
#endif
}
Output
Getting matrix rows on host:
Row: 0 -> 0:2
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 -> 4:2
Divide each matrix row by a sum of all elements in the row ...
Row: 0 -> 0:1
Row: 1 -> 0:0.25 1:0.5 2:0.25
Row: 2 -> 1:0.25 2:0.5 3:0.25
Row: 3 -> 2:0.25 3:0.5 4:0.25
Row: 4 -> 4:1
Getting matrix rows on CUDA device:
Row: 0 -> 0:2
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 -> 4:2
Divide each matrix row by a sum of all elements in the row ...
Row: 0 -> 0:1
Row: 1 -> 0:0.25 1:0.5 2:0.25
Row: 2 -> 1:0.25 2:0.5 3:0.25
Row: 3 -> 2:0.25 3:0.5 4:0.25
Row: 4 -> 4:1

◆ getColumnIndexes() [1/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
auto TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::getColumnIndexes

Getter of column indexes for nonconstant instances.

Returns
Reference to a vector with matrix elements column indexes.

◆ getColumnIndexes() [2/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
auto TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::getColumnIndexes

Getter of column indexes for constant instances.

Returns
Constant reference to a vector with matrix elements column indexes.

◆ getCompressedRowLengths()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Vector >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void getCompressedRowLengthsExample()
{
triangularMatrix.setElements( {
{ 0, 0, 1 },
{ 1, 0, 2 }, { 1, 1, 3 },
{ 2, 0, 4 }, { 2, 1, 5 }, { 2, 2, 6 },
{ 3, 0, 7 }, { 3, 1, 8 }, { 3, 2, 9 }, { 3, 3, 10 },
{ 4, 0, 11 }, { 4, 1, 12 }, { 4, 2, 13 }, { 4, 3, 14 }, { 4, 4, 15 } } );
std::cout << triangularMatrix << std::endl;
triangularMatrix.getCompressedRowLengths( rowLengths );
std::cout << "Compressed row lengths are: " << rowLengths << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Getting compressed row lengths on host: " << std::endl;
getCompressedRowLengthsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Getting compressed row lengths on CUDA device: " << std::endl;
getCompressedRowLengthsExample< TNL::Devices::Cuda >();
#endif
}
Output
Getting compressed row lengths on host:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:3
Row: 2 -> 0:4 1:5 2:6
Row: 3 -> 0:7 1:8 2:9 3:10
Row: 4 -> 0:11 1:12 2:13 3:14 4:15
Compressed row lengths are: [ 1, 2, 3, 4, 5 ]
Getting compressed row lengths on CUDA device:
Row: 0 -> 0:1
Row: 1 -> 0:2 1:3
Row: 2 -> 0:4 1:5 2:6
Row: 3 -> 0:7 1:8 2:9 3:10
Row: 4 -> 0:11 1:12 2:13 3:14 4:15
Compressed row lengths are: [ 1, 2, 3, 4, 5 ]

◆ getConstView()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
auto TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::getConstView

Returns a non-modifiable view of the sparse matrix.

See SparseMatrixView.

Returns
sparse matrix view.

◆ getElement()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
__cuda_callable__ auto TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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. SparseMatrix::getRow or SparseMatrix::forElements and SparseMatrix::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/SparseMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void getElements()
{
5, // number of matrix rows
5, // number of matrix columns
{ // matrix elements definition
{ 0, 0, 2.0 },
{ 1, 0, -1.0 }, { 1, 1, 2.0 }, { 1, 2, -1.0 },
{ 2, 1, -1.0 }, { 2, 2, 2.0 }, { 2, 3, -1.0 },
{ 3, 2, -1.0 }, { 3, 3, 2.0 }, { 3, 4, -1.0 },
{ 4, 4, 2.0 } } );
for( int i = 0; i < 5; i++ )
{
for( int j = 0; j < 5; 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 0 0 0 0
-1 2 -1 0 0
0 -1 2 -1 0
0 0 -1 2 -1
0 0 0 0 2
Get elements on CUDA device:
2 0 0 0 0
-1 2 -1 0 0
0 -1 2 -1 0
0 0 -1 2 -1
0 0 0 0 2

◆ getNonzeroElementsCount()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
Index TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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 , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
__cuda_callable__ Index TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::getPaddingIndex

Returns a padding index value.

Padding index is used for column indexes of padding zeros. Padding zeros are used in some sparse matrix formats for better data alignment in memory.

Returns
value of the padding index.

◆ getRow() [1/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
__cuda_callable__ auto TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
#include <TNL/Pointers/SharedPointer.h>
template< typename Device >
void getRowExample()
{
auto rowCapacities = { 1, 1, 1, 1, 1 }; // Variadic templates in SharedPointer
// constructor do not recognize initializer
// list so we give it a hint.
TNL::Pointers::SharedPointer< MatrixType > matrix( rowCapacities, 5 );
MatrixType* matrix_device = &matrix.template modifyData< Device >();
auto f = [=] __cuda_callable__ ( int rowIdx ) mutable {
auto row = matrix_device->getRow( rowIdx );
row.setElement( 0, rowIdx, 10 * ( rowIdx + 1 ) );
};
/***
* For the case when Device is CUDA device we need to synchronize smart
* pointers. To avoid this you may use SparseMatrixView. See
* SparseMatrixView::getRow example for details.
*/
TNL::Pointers::synchronizeSmartPointersOnDevice< Device >();
/***
* Set the matrix elements.
*/
TNL::Algorithms::parallelFor< Device >( 0, matrix->getRows(), f );
std::cout << *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:10
Row: 1 -> 1:20
Row: 2 -> 2:30
Row: 3 -> 3:40
Row: 4 -> 4:50
Getting matrix rows on CUDA device:
Row: 0 -> 0:10
Row: 1 -> 1:20
Row: 2 -> 2:30
Row: 3 -> 3:40
Row: 4 -> 4:50

See SparseMatrixRowView.

◆ getRow() [2/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
__cuda_callable__ auto TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
#include <TNL/Pointers/SharedPointer.h>
template< typename Device >
void getRowExample()
{
matrix->setElements(
{ { 0, 0, 1 },
{ 1, 0, 1 }, { 1, 1, 2 },
{ 2, 0, 1 }, { 2, 1, 2 }, { 2, 2, 3 },
{ 3, 0, 1 }, { 3, 1, 2 }, { 3, 2, 3 }, { 3, 3, 4 },
{ 4, 0, 1 }, { 4, 1, 2 }, { 4, 2, 3 }, { 4, 3, 4 }, { 4, 4, 5 } } );
/***
* 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( rowIdx );
};
/***
* For the case when Device is CUDA device we need to synchronize smart
* pointers. To avoid this you may use SparseMatrixView. See
* SparseMatrixView::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 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 trace is 15.
Getting matrix rows on CUDA device:
Matrix trace is 15.

See SparseMatrixRowView.

◆ getRowCapacities()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Vector >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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.

◆ getRowCapacity()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
__cuda_callable__ Index TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::getRowCapacity ( IndexType  row) const

Returns capacity of given matrix row.

Parameters
rowindex of matrix row.
Returns
number of matrix elements allocated for the row.

◆ getSegments() [1/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
auto TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::getSegments

Getter of segments for non-constant instances.

Segments are a structure for addressing the matrix elements columns and values. In fact, Segments represent the sparse matrix format.

Returns
Non-constant reference to segments.

◆ getSegments() [2/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
auto TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::getSegments

Getter of segments for constant instances.

Segments are a structure for addressing the matrix elements columns and values. In fact, Segments represent the sparse matrix format.

Returns
Constant reference to segments.

◆ getSerializationType()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
std::string TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::getSerializationType
static

Returns string with serialization type.

The string has a form Matrices::SparseMatrix< RealType, [any_device], IndexType, General/Symmetric, Format, [any_allocator] >.

Returns
String with the serialization type.
Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.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: SparseMatrix.hpp:193
Output
Get serialization type on CPU ...
Matrix type is: Matrices::SparseMatrix< double, CSR< [any_device], int, TNL::Algorithms::Segments::CSRScalarKernel<int, TNL::Devices::Host> >, [any_device], int, General, [any_allocator], [any_allocator] >Get serialization type on CUDA GPU ...
Matrix type is: Matrices::SparseMatrix< double, CSR< [any_device], int, TNL::Algorithms::Segments::CSRScalarKernel<int, TNL::Devices::Cuda> >, [any_device], int, General, [any_allocator], [any_allocator] >

◆ getSerializationTypeVirtual()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
std::string TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::getSerializationTypeVirtual
overridevirtual

Returns string with serialization type.

See SparseMatrix::getSerializationType.

Returns
String with the serialization type.
Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.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::SparseMatrix< double, CSR< [any_device], int, TNL::Algorithms::Segments::CSRScalarKernel<int, TNL::Devices::Host> >, [any_device], int, General, [any_allocator], [any_allocator] >Get serialization type on CUDA GPU ...
Matrix type is: Matrices::SparseMatrix< double, CSR< [any_device], int, TNL::Algorithms::Segments::CSRScalarKernel<int, TNL::Devices::Cuda> >, [any_device], int, General, [any_allocator], [any_allocator] >

Reimplemented from TNL::Object.

◆ getView()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
auto TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::getView

Returns a modifiable view of the sparse matrix.

See SparseMatrixView.

Returns
sparse matrix view.

◆ isBinary()

template<typename Real = double, typename Device = Devices::Host, typename Index = int, typename MatrixType = GeneralMatrix, template< typename Device_, typename Index_, typename IndexAllocator_ > class Segments = Algorithms::Segments::CSRDefault, typename ComputeReal = typename ChooseSparseMatrixComputeReal< Real, Index >::type, typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
static constexpr bool TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::isBinary ( )
inlinestaticconstexpr

Test of binary matrix type.

Returns
true if the matrix is stored as binary and false otherwise.

◆ isSymmetric()

template<typename Real = double, typename Device = Devices::Host, typename Index = int, typename MatrixType = GeneralMatrix, template< typename Device_, typename Index_, typename IndexAllocator_ > class Segments = Algorithms::Segments::CSRDefault, typename ComputeReal = typename ChooseSparseMatrixComputeReal< Real, Index >::type, typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
static constexpr bool TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::isSymmetric ( )
inlinestaticconstexpr

Test of symmetric matrix type.

Returns
true if the matrix is stored as symmetric and false otherwise.

◆ load() [1/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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 , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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 , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Matrix >
bool TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::operator!= ( const Matrix matrix) const

Comparison operator with another arbitrary matrix type.

Parameters
matrixis the right-hand side matrix.
Returns
true if the RHS matrix is equal, false otherwise.

◆ operator=() [1/4]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, typename MatrixType = GeneralMatrix, template< typename Device_, typename Index_, typename IndexAllocator_ > class Segments = Algorithms::Segments::CSRDefault, typename ComputeReal = typename ChooseSparseMatrixComputeReal< Real, Index >::type, 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_ >
SparseMatrix & TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::operator= ( const DenseMatrix< Real_, Device_, Index_, Organization, RealAllocator_ > &  matrix)

Assignment of dense matrix.

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

◆ operator=() [2/4]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, typename MatrixType = GeneralMatrix, template< typename Device_, typename Index_, typename IndexAllocator_ > class Segments = Algorithms::Segments::CSRDefault, typename ComputeReal = typename ChooseSparseMatrixComputeReal< Real, Index >::type, typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
template<typename RHSMatrix >
SparseMatrix & TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::operator= ( const RHSMatrix &  matrix)

Assignment of any matrix type other then this and dense.

Warning: Assignment of symmetric sparse matrix to general sparse matrix does not give correct result, currently. Only the diagonal and the lower part of the matrix is assigned.

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

◆ operator=() [3/4]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator > & TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::operator= ( const SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator > &  matrix)

Copy-assignment of exactly the same matrix type.

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

◆ operator=() [4/4]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator > & TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::operator= ( SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator > &&  matrix)
noexcept

Move-assignment of exactly the same matrix type.

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

◆ operator==()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Matrix >
bool TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::operator== ( const Matrix matrix) const

Comparison operator with another arbitrary matrix type.

Parameters
matrixis the right-hand side matrix.
Returns
true if the RHS matrix is equal, false otherwise.

◆ print()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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 , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::reduceAllRows ( Fetch &  fetch,
const 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 ) -> FetchValu { ... };
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/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
{
{ 0, 0, 1 },
{ 1, 1, 1 }, { 1, 2, 8 },
{ 2, 2, 1 }, { 2, 3, 9 },
{ 3, 3, 1 }, { 3, 4, 9 },
{ 4, 4, 1 } } );
/***
* Find largest element in each row.
*/
TNL::Containers::Vector< double, Device > rowMax( matrix.getRows() );
/***
* Prepare vector view and matrix 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 << "All rows reduction on host:" << std::endl;
reduceAllRows< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "All rows reduction on CUDA device:" << std::endl;
reduceAllRows< TNL::Devices::Cuda >();
#endif
}
void reduceAllRows(Fetch &fetch, const Reduce &reduce, Keep &keep, const FetchReal &identity)
Method for performing general reduction on all matrix rows.
Definition: SparseMatrix.hpp:625
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
__cuda_callable__ T abs(const T &n)
This function returns absolute value of given number n.
Definition: Math.h:91
constexpr ResultType max(const T1 &a, const T2 &b)
This function returns maximum of two numbers.
Definition: Math.h:65
Output
All rows reduction on host:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 1:1 2:8
Row: 2 -> 2:1 3:9
Row: 3 -> 3:1 4:9
Row: 4 -> 4:1
Max. elements in rows are: [ 1, 8, 9, 9, 1 ]
All rows reduction on CUDA device:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 1:1 2:8
Row: 2 -> 2:1 3:9
Row: 3 -> 3:1 4:9
Row: 4 -> 4:1
Max. elements in rows are: [ 1, 8, 9, 9, 1 ]

◆ reduceAllRows() [2/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::reduceAllRows ( Fetch &  fetch,
const 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/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
{
{ 0, 0, 1 },
{ 1, 1, 1 }, { 1, 2, 8 },
{ 2, 2, 1 }, { 2, 3, 9 },
{ 3, 3, 1 }, { 3, 4, 9 },
{ 4, 4, 1 } } );
/***
* Find largest element in each row.
*/
TNL::Containers::Vector< double, Device > rowMax( matrix.getRows() );
/***
* Prepare vector view and matrix 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 << "All rows reduction on host:" << std::endl;
reduceAllRows< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "All rows reduction on CUDA device:" << std::endl;
reduceAllRows< TNL::Devices::Cuda >();
#endif
}
Output
All rows reduction on host:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 1:1 2:8
Row: 2 -> 2:1 3:9
Row: 3 -> 3:1 4:9
Row: 4 -> 4:1
Max. elements in rows are: [ 1, 8, 9, 9, 1 ]
All rows reduction on CUDA device:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 1:1 2:8
Row: 2 -> 2:1 3:9
Row: 3 -> 3:1 4:9
Row: 4 -> 4:1
Max. elements in rows are: [ 1, 8, 9, 9, 1 ]

◆ reduceRows() [1/2]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, typename MatrixType = GeneralMatrix, template< typename Device_, typename Index_, typename IndexAllocator_ > class Segments = Algorithms::Segments::CSRDefault, typename ComputeReal = typename ChooseSparseMatrixComputeReal< Real, Index >::type, typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::reduceRows ( IndexType  begin,
IndexType  end,
Fetch &  fetch,
const 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/SparseMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void reduceRows()
{
{ 0, 0, 1 },
{ 1, 1, 1 }, { 1, 2, 8 },
{ 2, 2, 1 }, { 2, 3, 9 },
{ 3, 3, 1 }, { 3, 4, 9 },
{ 4, 4, 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.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, const Reduce &reduce, Keep &keep, const FetchReal &identity)
Method for performing general reduction on matrix rows.
Output
Rows reduction on host:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 1:1 2:8
Row: 2 -> 2:1 3:9
Row: 3 -> 3:1 4:9
Row: 4 -> 4:1
Max. elements in rows are: [ 1, 8, 9, 9, 1 ]
Rows reduction on CUDA device:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 1:1 2:8
Row: 2 -> 2:1 3:9
Row: 3 -> 3:1 4:9
Row: 4 -> 4:1
Max. elements in rows are: [ 1, 8, 9, 9, 1 ]

◆ reduceRows() [2/2]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, typename MatrixType = GeneralMatrix, template< typename Device_, typename Index_, typename IndexAllocator_ > class Segments = Algorithms::Segments::CSRDefault, typename ComputeReal = typename ChooseSparseMatrixComputeReal< Real, Index >::type, typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::reduceRows ( IndexType  begin,
IndexType  end,
Fetch &  fetch,
const 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_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/SparseMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void reduceRows()
{
{ 0, 0, 1 },
{ 1, 1, 1 }, { 1, 2, 8 },
{ 2, 2, 1 }, { 2, 3, 9 },
{ 3, 3, 1 }, { 3, 4, 9 },
{ 4, 4, 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.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 -> 1:1 2:8
Row: 2 -> 2:1 3:9
Row: 3 -> 3:1 4:9
Row: 4 -> 4:1
Max. elements in rows are: [ 1, 8, 9, 9, 1 ]
Rows reduction on CUDA device:
The matrix reads as:
Row: 0 -> 0:1
Row: 1 -> 1:1 2:8
Row: 2 -> 2:1 3:9
Row: 3 -> 3:1 4:9
Row: 4 -> 4:1
Max. elements in rows are: [ 1, 8, 9, 9, 1 ]

◆ save() [1/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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 , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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 , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::sequentialForAllRows ( Function &  function)

This method calls sequentialForRows for all matrix rows.

See SparseMatrix::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 , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::sequentialForAllRows ( Function &  function) const

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

See SparseMatrix::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 , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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__ ( RowView& row ) { ... };

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

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 , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Function >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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 = [] __cuda_callable__ ( RowView& row ) { ... };

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

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.

◆ setColumnsWithoutReset()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::setColumnsWithoutReset ( IndexType  columns)
virtual

Set number of columns of this matrix.

Unlike setDimensions, the storage is not reset in this operation. It is the user's responsibility to update the column indices stored in the matrix to be consistent with the new number of columns.

Parameters
columnsis the number of matrix columns.

◆ setDimensions()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::setDimensions ( IndexType  rows,
IndexType  columns 
)
overridevirtual

Set number of rows and columns of this matrix.

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

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

◆ setElement()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
__cuda_callable__ void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, 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. SparseMatrix::getRow or SparseMatrix::forElements and SparseMatrix::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/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
#include <TNL/Pointers/SharedPointer.h>
#include <TNL/Pointers/SmartPointersRegister.h>
template< typename Device >
{
auto rowCapacities = { 1, 1, 1, 1, 1 };
TNL::Pointers::SharedPointer< MatrixType > matrix( rowCapacities, 5 );
/***
* Calling the method setElements from host (CPU).
*/
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;
/***
* This lambda function will run on the native device of the matrix which can be CPU or GPU.
*/
MatrixType* matrix_device = &matrix.template modifyData< Device >();
auto f = [=] __cuda_callable__ ( int i ) mutable {
matrix_device->setElement( i, i, -i );
};
/***
* For the case when Device is CUDA device we need to synchronize smart
* pointers. To avoid this you may use SparseMatrixView. See
* SparseMatrixView::getRow example for details.
*/
TNL::Pointers::synchronizeSmartPointersOnDevice< Device >();
TNL::Algorithms::parallelFor< Device >( 0, 5, 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
}
void setElements(const std::initializer_list< std::tuple< IndexType, IndexType, RealType > > &data)
This method sets the sparse matrix elements from initializer list.
Definition: SparseMatrix.hpp:330
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 ->
Row: 1 -> 1:-1
Row: 2 -> 2:-2
Row: 3 -> 3:-3
Row: 4 -> 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 ->
Row: 1 -> 1:-1
Row: 2 -> 2:-2
Row: 3 -> 3:-3
Row: 4 -> 4:-4

◆ setElements() [1/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::setElements ( const std::initializer_list< std::tuple< IndexType, IndexType, RealType > > &  data)

This method sets the sparse matrix elements from initializer list.

The number of matrix rows and columns must be set already. The matrix elements values are given as a list data of triples: { { row1, column1, value1 }, { row2, column2, value2 }, ... }.

Parameters
datais a initializer list of initializer lists representing list of matrix rows.
Example
#include <iostream>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void setElementsExample()
{
TNL::Matrices::SparseMatrix< double, Device > matrix ( 5, 5 ); // matrix dimensions
matrix.setElements( { // matrix elements definition
{ 0, 0, 2.0 },
{ 1, 0, -1.0 }, { 1, 1, 2.0 }, { 1, 2, -1.0 },
{ 2, 1, -1.0 }, { 2, 2, 2.0 }, { 2, 3, -1.0 },
{ 3, 2, -1.0 }, { 3, 3, 2.0 }, { 3, 4, -1.0 },
{ 4, 4, 2.0 } } );
std::cout << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Setting matrix elements on host: " << std::endl;
setElementsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Setting matrix elements on CUDA device: " << std::endl;
setElementsExample< TNL::Devices::Cuda >();
#endif
}
Output
Setting matrix elements on host:
Row: 0 -> 0:2
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 -> 4:2
Setting matrix elements on CUDA device:
Row: 0 -> 0:2
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 -> 4:2

◆ setElements() [2/2]

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename MapIndex , typename MapValue >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::setElements ( const std::map< std::pair< MapIndex, MapIndex >, MapValue > &  map)

This method sets the sparse matrix elements from std::map.

The matrix elements values are given as a map data where keys are std::pair of matrix coordinates ( {row, column} ) and value is the matrix element value.

Template Parameters
MapIndexis a type for indexing rows and columns.
MapValueis a type for matrix elements values in the map.
Parameters
mapis std::map containing matrix elements.
Example
#include <iostream>
#include <map>
#include <utility>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void setElementsExample()
{
map.insert( std::make_pair( std::make_pair( 0, 0 ), 2.0 ) );
map.insert( std::make_pair( std::make_pair( 1, 0 ), -1.0 ) );
map.insert( std::make_pair( std::make_pair( 1, 1 ), 2.0 ) );
map.insert( std::make_pair( std::make_pair( 1, 2 ), -1.0 ) );
map.insert( std::make_pair( std::make_pair( 2, 1 ), -1.0 ) );
map.insert( std::make_pair( std::make_pair( 2, 2 ), 2.0 ) );
map.insert( std::make_pair( std::make_pair( 2, 3 ), -1.0 ) );
map.insert( std::make_pair( std::make_pair( 3, 2 ), -1.0 ) );
map.insert( std::make_pair( std::make_pair( 3, 3 ), 2.0 ) );
map.insert( std::make_pair( std::make_pair( 3, 4 ), -1.0 ) );
map.insert( std::make_pair( std::make_pair( 4, 4 ), 2.0 ) );
matrix.setElements( map );
std::cout << "General sparse matrix: " << std::endl << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrices on CPU ... " << std::endl;
setElementsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrices on CUDA GPU ... " << std::endl;
setElementsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrices on CPU ...
General sparse matrix:
Row: 0 -> 0:2
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 -> 4:2
Creating matrices on CUDA GPU ...
General sparse matrix:
Row: 0 -> 0:2
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 -> 4:2

◆ setLike()

template<typename Real = double, typename Device = Devices::Host, typename Index = int, typename MatrixType = GeneralMatrix, template< typename Device_, typename Index_, typename IndexAllocator_ > class Segments = Algorithms::Segments::CSRDefault, typename ComputeReal = typename ChooseSparseMatrixComputeReal< Real, Index >::type, typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >, typename IndexAllocator = typename Allocators::Default< Device >::template Allocator< Index >>
template<typename Matrix >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::setLike ( const Matrix matrix)

Set the number of matrix rows and columns by the given matrix.

Template Parameters
Matrixis matrix type. This can be any matrix having methods getRows and getColumns.
Parameters
matrixin the input matrix dimensions of which are to be adopted.

◆ setRowCapacities()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename RowsCapacitiesVector >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::setRowCapacities ( const RowsCapacitiesVector &  rowCapacities)

Allocates memory for non-zero matrix elements.

The size of the input vector must be equal to the number of matrix rows. The number of allocated matrix elements for each matrix row depends on the sparse matrix format. Some formats may allocate more elements than required.

Template Parameters
RowsCapacitiesVectoris a type of vector/array used for row capacities setting.
Parameters
rowCapacitiesis a vector telling the number of required non-zero matrix elements in each row.
Example
#include <iostream>
#include <TNL/Containers/Vector.h>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void setRowCapacitiesExample()
{
TNL::Containers::Vector< int, Device > rowCapacities{ 1, 2, 3, 4, 5 };
matrix.setRowCapacities( rowCapacities );
for( int row = 0; row < 5; row++ )
for( int column = 0; column <= row; column++ )
matrix.setElement( row, column, row - column + 1 );
std::cout << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrix on CPU ... " << std::endl;
setRowCapacitiesExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA GPU ... " << std::endl;
setRowCapacitiesExample< TNL::Devices::Cuda >();
#endif
}
__cuda_callable__ void setElement(IndexType i, ValueType value)
Sets the value of the i-th element to v.
Definition: Array.hpp:356
Output
Creating matrix on CPU ...
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 0:4 1:3 2:2 3:1
Row: 4 -> 0:5 1:4 2:3 3:2 4:1
Creating matrix on CUDA GPU ...
Row: 0 -> 0:1
Row: 1 -> 0:2 1:1
Row: 2 -> 0:3 1:2 2:1
Row: 3 -> 0:4 1:3 2:2 3:1
Row: 4 -> 0:5 1:4 2:3 3:2 4:1

◆ vectorProduct()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename InVector , typename OutVector >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::vectorProduct ( const InVector &  inVector,
OutVector &  outVector,
ComputeRealType  matrixMultiplicator = 1.0,
ComputeRealType  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: