Template Numerical Library version\ main:94209208
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Protected Attributes | List of all members
TNL::Matrices::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 ColumnIndexesVectorType = Containers::Vector< Index, Device, Index, IndexAllocator >
 Type of vector holding values of column indexes.
 
using ConstViewType
 Matrix view type for constant instances.
 
using IndexAllocatorType = IndexAllocator
 The allocator for matrix elements column indexes.
 
using MatrixType = MatrixType_
 The type of matrix - general, symmetric or binary.
 
using RealAllocatorType = RealAllocator
 The allocator for matrix elements values.
 
using RowCapacitiesVectorType = Containers::Vector< Index, Device, Index, IndexAllocator >
 Type of vector holding values of row capacities.
 
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.
 
template<typename _Real = Real, typename _Device = Device, typename _Index = Index, typename _MatrixType = MatrixType, template< typename, typename, typename > class _Segments = SegmentsTemplate, 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 = Containers::Vector< Real, Device, Index, RealAllocator >
 Type of vector holding values of matrix elements.
 
using ViewType = SparseMatrixView< Real, Device, Index, MatrixType, SegmentsViewTemplate, ComputeReal >
 Type of related matrix view.
 
- Public Types inherited from TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >
using ColumnIndexesViewType
 
using ComputeRealType = ComputeReal
 
using ConstRowView = typename RowView::ConstRowView
 Type for accessing constant matrix rows.
 
using DefaultSegmentsReductionKernel = typename Algorithms::SegmentsReductionKernels::DefaultKernel< SegmentsView >::type
 Type of the kernel used for parallel reductions on segments.
 
using DeviceType = Device
 The device where the matrix is allocated.
 
using IndexType = Index
 The type used for matrix elements indexing.
 
using RealType = typename Base::RealType
 The type of matrix elements.
 
using RowView
 Type for accessing matrix rows.
 
using SegmentsViewType = SegmentsView
 Type of segments view used by this matrix. It represents the sparse matrix format.
 
- Public Types inherited from TNL::Matrices::MatrixBase< Real, Device, Index, MatrixType, SegmentsView::getOrganization() >
using ConstValuesViewType
 Type of constant vector view holding values of matrix elements.
 
using DeviceType
 The device where the matrix is allocated.
 
using IndexType
 The type used for matrix elements indexing.
 
using RealType
 The type of matrix elements.
 
using RowCapacitiesType
 
using ValuesViewType
 Type of vector view holding values of matrix elements.
 

Public Member Functions

 SparseMatrix (const RealAllocatorType &realAllocator=RealAllocatorType(), const IndexAllocatorType &indexAllocator=IndexAllocatorType())
 Constructor only with values and column indexes allocators.
 
template<typename RowCapacitiesVector , std::enable_if_t< TNL::IsArrayType< RowCapacitiesVector >::value, int > = 0>
 SparseMatrix (const RowCapacitiesVector &rowCapacities, Index columns, const RealAllocatorType &realAllocator=RealAllocatorType(), const IndexAllocatorType &indexAllocator=IndexAllocatorType())
 Constructor with matrix rows capacities given as a vector and number of columns.
 
 SparseMatrix (const SparseMatrix &matrix)
 Copy constructor.
 
template<typename ListIndex >
 SparseMatrix (const std::initializer_list< ListIndex > &rowCapacities, Index columns, const RealAllocatorType &realAllocator=RealAllocatorType(), const IndexAllocatorType &indexAllocator=IndexAllocatorType())
 Constructor with matrix rows capacities and number of columns.
 
 SparseMatrix (Index rows, Index columns, const std::initializer_list< std::tuple< Index, Index, Real > > &data, SymmetricMatrixEncoding encoding=SymmetricMatrixEncoding::LowerPart, const RealAllocatorType &realAllocator=RealAllocatorType(), const IndexAllocatorType &indexAllocator=IndexAllocatorType())
 Constructor with matrix dimensions and data in initializer list.
 
template<typename MapIndex , typename MapValue >
 SparseMatrix (Index rows, Index columns, const std::map< std::pair< MapIndex, MapIndex >, MapValue > &map, SymmetricMatrixEncoding encoding=SymmetricMatrixEncoding::LowerPart, const RealAllocatorType &realAllocator=RealAllocatorType(), const IndexAllocatorType &indexAllocator=IndexAllocatorType())
 Constructor with matrix dimensions and data in std::map.
 
template<typename Index_t , std::enable_if_t< std::is_integral_v< Index_t >, int > = 0>
 SparseMatrix (Index_t rows, Index_t columns, const RealAllocatorType &realAllocator=RealAllocatorType(), const IndexAllocatorType &indexAllocator=IndexAllocatorType())
 Constructor with matrix dimensions.
 
 SparseMatrix (SparseMatrix &&) noexcept=default
 Move constructor.
 
ConstViewType getConstView () const
 Returns a non-modifiable view of the sparse matrix.
 
template<typename Real2 , typename Index2 , template< typename, typename, typename > class Segments2>
void getTransposition (const SparseMatrix< Real2, Device, Index2, MatrixType, Segments2 > &matrix, const ComputeReal &matrixMultiplicator=1.0)
 Computes transposition of the matrix.
 
ViewType getView ()
 Returns a modifiable view of the sparse matrix.
 
void load (const String &fileName)
 Method for loading the matrix from the file with given filename.
 
void load (File &file) override
 Method for loading the matrix from a file.
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization, typename RealAllocator_ >
SparseMatrixoperator= (const DenseMatrix< Real_, Device_, Index_, Organization, RealAllocator_ > &matrix)
 Assignment of dense matrix.
 
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.
 
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.
 
SparseMatrixoperator= (SparseMatrix &&matrix) noexcept(false)
 Move-assignment of exactly the same matrix type.
 
void reset ()
 Resets the matrix to zero dimensions.
 
void save (const String &fileName) const
 Method for saving the matrix to the file with given filename.
 
void save (File &file) const override
 Method for saving the matrix to a file.
 
virtual void setColumnsWithoutReset (Index columns)
 Set number of columns of this matrix.
 
void setDimensions (Index rows, Index columns)
 Set number of rows and columns of this matrix.
 
void setElements (const std::initializer_list< std::tuple< Index, Index, Real > > &data, SymmetricMatrixEncoding encoding=SymmetricMatrixEncoding::LowerPart)
 This method sets the sparse matrix elements from initializer list.
 
template<typename MapIndex , typename MapValue >
void setElements (const std::map< std::pair< MapIndex, MapIndex >, MapValue > &map, SymmetricMatrixEncoding encoding=SymmetricMatrixEncoding::LowerPart)
 This method sets the sparse matrix elements from std::map.
 
template<typename Matrix >
void setLike (const Matrix &matrix)
 Set the number of matrix rows and columns by the given matrix.
 
template<typename Matrix_ >
void setLike (const Matrix_ &matrix)
 
template<typename RowCapacitiesVector >
void setRowCapacities (const RowCapacitiesVector &rowCapacities)
 Allocates memory for non-zero matrix elements.
 
- Public Member Functions inherited from TNL::Object
virtual ~Object ()=default
 Destructor.
 
virtual std::string getSerializationTypeVirtual () const
 
void load (const String &fileName)
 Method for restoring the object from a file.
 
void save (const String &fileName) const
 Method for saving the object to a file as a binary data.
 
- Public Member Functions inherited from TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >
__cuda_callable__ SparseMatrixBase ()=default
 Constructor with no parameters.
 
__cuda_callable__ SparseMatrixBase (const SparseMatrixBase &matrix)=default
 Copy constructor.
 
__cuda_callable__ SparseMatrixBase (IndexType rows, IndexType columns, typename Base::ValuesViewType values, ColumnIndexesViewType columnIndexes, SegmentsViewType segments)
 Constructor with all necessary data and views.
 
__cuda_callable__ SparseMatrixBase (SparseMatrixBase &&matrix) noexcept=default
 Move constructor.
 
__cuda_callable__ void addElement (IndexType row, IndexType column, const RealType &value, const RealType &thisElementMultiplicator=1.0)
 Add element at given row and column to given value.
 
__cuda_callable__ IndexType findElement (IndexType row, IndexType column) const
 Finds element in the matrix and returns its position in the arrays with values and columnIndexes.
 
template<typename Function >
void forAllElements (Function &&function)
 This method calls forElements for all matrix rows.
 
template<typename Function >
void forAllElements (Function &&function) const
 This method calls forElements for all matrix rows (for constant instances).
 
template<typename Function >
void forAllRows (Function &&function)
 Method for parallel iteration over all matrix rows.
 
template<typename Function >
void forAllRows (Function &&function) const
 Method for parallel iteration over all matrix rows for constant instances.
 
template<typename Function >
void forElements (IndexType begin, IndexType end, Function &&function)
 Method for iteration over all matrix rows for non-constant instances.
 
template<typename Function >
void forElements (IndexType begin, IndexType end, Function &&function) const
 Method for iteration over all matrix rows for constant instances.
 
template<typename Function >
void forRows (IndexType begin, IndexType end, Function &&function)
 Method for parallel iteration over matrix rows from interval [begin, end).
 
template<typename Function >
void forRows (IndexType begin, IndexType end, Function &&function) const
 Method for parallel iteration over matrix rows from interval [begin, end) for constant instances.
 
ColumnIndexesViewTypegetColumnIndexes ()
 Getter of column indexes for nonconstant instances.
 
const ColumnIndexesViewTypegetColumnIndexes () const
 Getter of column indexes for constant instances.
 
template<typename Vector >
void getCompressedRowLengths (Vector &rowLengths) const
 Computes number of non-zeros in each row.
 
__cuda_callable__ RealType getElement (IndexType row, IndexType column) const
 Returns value of matrix element at position given by its row and column index.
 
IndexType getNonzeroElementsCount () const override
 Returns number of non-zero matrix elements.
 
__cuda_callable__ RowView getRow (IndexType rowIdx)
 Non-constant getter of simple structure for accessing given matrix row.
 
__cuda_callable__ ConstRowView getRow (IndexType rowIdx) const
 Constant getter of simple structure for accessing given matrix row.
 
template<typename Vector >
void getRowCapacities (Vector &rowCapacities) const
 Compute capacities of all rows.
 
__cuda_callable__ IndexType getRowCapacity (IndexType row) const
 Returns capacity of given matrix row.
 
SegmentsViewTypegetSegments ()
 Getter of segments for non-constant instances.
 
const SegmentsViewTypegetSegments () const
 Getter of segments for constant instances.
 
template<typename Matrix >
bool operator!= (const Matrix &matrix) const
 Comparison operator with another arbitrary matrix type.
 
SparseMatrixBaseoperator= (const SparseMatrixBase &)=delete
 Copy-assignment operator.
 
SparseMatrixBaseoperator= (SparseMatrixBase &&)=delete
 Move-assignment operator.
 
template<typename Matrix >
bool operator== (const Matrix &matrix) const
 Comparison operator with another arbitrary matrix type.
 
void print (std::ostream &str) const
 Method for printing the matrix to output stream.
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchValue , typename SegmentsReductionKernel = DefaultSegmentsReductionKernel>
void reduceAllRows (Fetch &fetch, const Reduce &reduce, Keep &keep, const FetchValue &identity, const SegmentsReductionKernel &kernel=SegmentsReductionKernel{}) const
 Method for performing general reduction on all matrix rows for constant instances.
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchValue , typename SegmentsReductionKernel = DefaultSegmentsReductionKernel>
void reduceRows (IndexType begin, IndexType end, Fetch &fetch, const Reduce &reduce, Keep &keep, const FetchValue &identity, const SegmentsReductionKernel &kernel=SegmentsReductionKernel{}) const
 Method for performing general reduction on matrix rows for constant instances.
 
template<typename Function >
void sequentialForAllRows (Function &&function)
 This method calls sequentialForRows for all matrix rows.
 
template<typename Function >
void sequentialForAllRows (Function &&function) const
 This method calls sequentialForRows for all matrix rows (for constant instances).
 
template<typename Function >
void sequentialForRows (IndexType begin, IndexType end, Function &&function)
 Method for sequential iteration over all matrix rows for non-constant instances.
 
template<typename Function >
void sequentialForRows (IndexType begin, IndexType end, Function &&function) const
 Method for sequential iteration over all matrix rows for constant instances.
 
__cuda_callable__ void setElement (IndexType row, IndexType column, const RealType &value)
 Sets element at given row and column to given value.
 
void sortColumnIndexes ()
 Sort matrix elements in each row by column indexes in ascending order.
 
template<typename InVector , typename OutVector >
void transposedVectorProduct (const InVector &inVector, OutVector &outVector, ComputeReal matrixMultiplicator=1.0, ComputeReal outVectorMultiplicator=0.0, Index begin=0, Index end=0) const
 Computes product of transposed matrix and vector.
 
template<typename InVector , typename OutVector , typename SegmentsReductionKernel = DefaultSegmentsReductionKernel>
void vectorProduct (const InVector &inVector, OutVector &outVector, ComputeRealType matrixMultiplicator=1.0, ComputeRealType outVectorMultiplicator=0.0, IndexType begin=0, IndexType end=0, const SegmentsReductionKernel &kernel=SegmentsReductionKernel{}) const
 Computes product of matrix and vector.
 
template<typename InVector , typename OutVector , typename SegmentsReductionKernel , typename... , std::enable_if_t< ! std::is_convertible_v< SegmentsReductionKernel, ComputeRealType >, bool > = true>
void vectorProduct (const InVector &inVector, OutVector &outVector, const SegmentsReductionKernel &kernel) const
 
- Public Member Functions inherited from TNL::Matrices::MatrixBase< Real, Device, Index, MatrixType, SegmentsView::getOrganization() >
__cuda_callable__ MatrixBase ()=default
 Basic constructor with no parameters.
 
__cuda_callable__ MatrixBase (const MatrixBase &view)=default
 Shallow copy constructor.
 
__cuda_callable__ MatrixBase (IndexType rows, IndexType columns, ValuesViewType values)
 Constructor with matrix dimensions and matrix elements values.
 
__cuda_callable__ MatrixBase (MatrixBase &&view) noexcept=default
 Move constructor.
 
IndexType getAllocatedElementsCount () const
 Tells the number of allocated matrix elements.
 
__cuda_callable__ IndexType getColumns () const
 Returns number of matrix columns.
 
__cuda_callable__ IndexType getRows () const
 Returns number of matrix rows.
 
__cuda_callable__ ValuesViewTypegetValues ()
 Returns a reference to a vector with the matrix elements values.
 
__cuda_callable__ const ValuesViewTypegetValues () const
 Returns a constant reference to a vector with the matrix elements values.
 
bool operator!= (const Matrix &matrix) const
 Comparison operator with another arbitrary matrix view type.
 
bool operator!= (const MatrixT &matrix) const
 
__cuda_callable__ MatrixBaseoperator= (const MatrixBase &)=delete
 Copy-assignment operator.
 
__cuda_callable__ MatrixBaseoperator= (MatrixBase &&)=delete
 Move-assignment operator.
 
bool operator== (const Matrix &matrix) const
 Comparison operator with another arbitrary matrix view type.
 
bool operator== (const MatrixT &matrix) const
 

Static Public Member Functions

static std::string getSerializationType ()
 Returns string with serialization type.
 
- Static Public Member Functions inherited from TNL::Object
static std::string getSerializationType ()
 Static serialization type getter.
 
- Static Public Member Functions inherited from TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >
static std::string getSerializationType ()
 Returns string with serialization type.
 
- Static Public Member Functions inherited from TNL::Matrices::MatrixBase< Real, Device, Index, MatrixType, SegmentsView::getOrganization() >
static constexpr ElementsOrganization getOrganization ()
 Matrix elements organization getter.
 
static constexpr bool isBinary ()
 Test of binary matrix type.
 
static constexpr bool isMatrix ()
 Test of matrix type.
 
static constexpr bool isSymmetric ()
 Test of symmetric matrix type.
 

Protected Attributes

ColumnIndexesVectorType columnIndexes
 Vector containing the column indexes of allocated elements.
 
SegmentsType segments
 Instance of the segments used for indexing in the sparse matrix.
 
ValuesVectorType values
 Vector containing the allocated matrix elements.
 
- Protected Attributes inherited from TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >
ColumnIndexesViewType columnIndexes
 
SegmentsViewType segments
 
- Protected Attributes inherited from TNL::Matrices::MatrixBase< Real, Device, Index, MatrixType, SegmentsView::getOrganization() >
IndexType columns
 
IndexType rows
 
ValuesViewType values
 

Additional Inherited Members

- Protected Member Functions inherited from TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >
__cuda_callable__ void bind (IndexType rows, IndexType columns, typename Base::ValuesViewType values, ColumnIndexesViewType columnIndexes, SegmentsViewType segments)
 Re-initializes the internal attributes of the base class.
 
- Protected Member Functions inherited from TNL::Matrices::MatrixBase< Real, Device, Index, MatrixType, SegmentsView::getOrganization() >
__cuda_callable__ void bind (IndexType rows, IndexType columns, ValuesViewType values)
 Re-initializes the internal attributes of the base class.
 

Detailed Description

template<typename Real = double, typename Device = Devices::Host, typename Index = int, typename MatrixType_ = GeneralMatrix, template< typename Device_, typename Index_, typename IndexAllocator_ > class Segments = Algorithms::Segments::CSR, 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::CSR, 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
Initial value:
SparseMatrixView< std::add_const_t< Real >, Device, Index, MatrixType, SegmentsViewTemplate, ComputeReal >
MatrixType_ MatrixType
The type of matrix - general, symmetric or binary.
Definition SparseMatrix.h:84
typename SegmentsType::template ViewTemplate< Device_, Index > SegmentsViewTemplate
Templated view type of segments, i.e. sparse matrix format.
Definition SparseMatrix.h:101

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::CSR, 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, ComputeReal >

Type of related matrix view.

See SparseMatrixView.

Constructor & Destructor Documentation

◆ SparseMatrix() [1/6]

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

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_v< Index_t >, 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() [3/6]

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,
Index 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()
{
TNL::Matrices::SparseMatrix< double, Device > matrix{ { 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
}
__cuda_callable__ void setElement(IndexType row, IndexType column, const RealType &value)
Sets element at given row and column to given value.
Definition SparseMatrixBase.hpp:154
Implementation of sparse matrix, i.e. matrix storing only non-zero elements.
Definition SparseMatrix.h:57
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() [4/6]

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,
Index 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 };
TNL::Matrices::SparseMatrix< double, Device > matrix{ 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:36
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() [5/6]

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 ( Index rows,
Index columns,
const std::initializer_list< std::tuple< Index, Index, Real > > & data,
SymmetricMatrixEncoding encoding = SymmetricMatrixEncoding::LowerPart,
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.
encodingdefines encoding for sparse symmetric matrices - see TNL::Matrices::SymmetricMatrixEncoding.
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()
{
// number of matrix rows
5,
// number of matrix columns
5,
// matrix elements definition
{
// clang-format off
{ 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 }
// clang-format on
} );
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() [6/6]

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 ( Index rows,
Index columns,
const std::map< std::pair< MapIndex, MapIndex >, MapValue > & map,
SymmetricMatrixEncoding encoding = SymmetricMatrixEncoding::LowerPart,
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.
encodingdefines encoding for sparse symmetric matrices - see TNL::Matrices::SymmetricMatrixEncoding.
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

◆ 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 ( ) const

Returns a non-modifiable view of the sparse matrix.

See SparseMatrixView.

Returns
sparse matrix view.

◆ getSerializationType()

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::CSR, 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 >>
std::string TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::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 SparseMatrixBase.hpp:43
Output
Get serialization type on CPU ...
Matrix type is: Matrices::SparseMatrix< double, CSR< int >, [any_device], int, General, [any_allocator], [any_allocator] >Get serialization type on CUDA GPU ...
Matrix type is: Matrices::SparseMatrix< double, CSR< int >, [any_device], int, General, [any_allocator], [any_allocator] >

◆ getTransposition()

template<typename Real , typename Device , typename Index , typename MatrixType , template< typename, typename, typename > class Segments, typename ComputeReal , typename RealAllocator , typename IndexAllocator >
template<typename Real2 , typename Index2 , template< typename, typename, typename > class Segments2>
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::getTransposition ( const SparseMatrix< Real2, Device, Index2, MatrixType, Segments2 > & matrix,
const ComputeReal & matrixMultiplicator = 1.0 )

Computes transposition of the matrix.

Template Parameters
Real2is the real type of the input matrix.
Index2is the index type of the input matrix.
Segments2is the type of the segments of the input matrix.
Parameters
matrixis the input matrix.
matrixMultiplicatoris the factor by which the matrix is multiplied.

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

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

◆ 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::CSR, 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::CSR, 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)

Move-assignment of exactly the same matrix type.

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

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

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

Set number of rows and columns of this matrix.

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

◆ 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< Index, Index, Real > > & data,
SymmetricMatrixEncoding encoding = SymmetricMatrixEncoding::LowerPart )

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.
encodingdefines encoding for sparse symmetric matrices - see TNL::Matrices::SymmetricMatrixEncoding.
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 elements definition
matrix.setElements( {
// clang-format off
{ 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 }
// clang-format on
} );
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,
SymmetricMatrixEncoding encoding = SymmetricMatrixEncoding::LowerPart )

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.
encodingdefines encoding for sparse symmetric matrices - see TNL::Matrices::SymmetricMatrixEncoding.
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::CSR, 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 RowCapacitiesVector >
void TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType, Segments, ComputeReal, RealAllocator, IndexAllocator >::setRowCapacities ( const RowCapacitiesVector & 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
RowCapacitiesVectoris 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:357
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

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