|
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.
|
using | ColumnIndexesViewType |
using | ComputeRealType |
using | ConstRowView |
| Type for accessing constant matrix rows.
|
using | DefaultSegmentsReductionKernel |
| Type of the kernel used for parallel reductions on segments.
|
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 | RowView |
| Type for accessing matrix rows.
|
using | SegmentsViewType |
| Type of segments view used by this matrix. It represents the sparse matrix format.
|
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.
|
|
| 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.
|
template<typename Real_, typename Device_, typename Index_, ElementsOrganization Organization, typename RealAllocator_> |
SparseMatrix & | operator= (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> |
SparseMatrix & | operator= (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) |
SparseMatrix & | operator= (const SparseMatrix &matrix) |
| Copy-assignment of exactly the same matrix type.
|
SparseMatrix & | operator= (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 | 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 | setDimensions (Index rows, Index columns, SegmentsType segments) |
| Set number of rows and columns of this matrix and performs full allocation according to the given segments.
|
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.
|
__cuda_callable__ | SparseMatrixBase ()=default |
| Constructor with no parameters.
|
__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.
|
void | forAllElements (Function &&function) const |
| This method calls forElements for all matrix rows (for constant instances).
|
void | forAllRows (Function &&function) |
| Method for parallel iteration over all matrix rows.
|
void | forElements (IndexType begin, IndexType end, Function &&function) const |
| Method for iteration over all matrix rows for constant instances.
|
void | forRows (IndexType begin, IndexType end, Function &&function) |
| Method for parallel iteration over matrix rows from interval [begin, end).
|
const ColumnIndexesViewType & | getColumnIndexes () const |
| Getter of column indexes for constant instances.
|
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 |
| Returns number of non-zero matrix elements.
|
__cuda_callable__ ConstRowView | getRow (IndexType rowIdx) const |
| Constant getter of simple structure for accessing given matrix row.
|
void | getRowCapacities (Vector &rowCapacities) const |
| Compute capacities of all rows.
|
__cuda_callable__ IndexType | getRowCapacity (IndexType row) const |
| Returns capacity of given matrix row.
|
SegmentsViewType & | getSegments () |
| Getter of segments for non-constant instances.
|
bool | operator!= (const Matrix &matrix) const |
| Comparison operator with another arbitrary matrix type.
|
SparseMatrixBase & | operator= (const SparseMatrixBase &)=delete |
| Copy-assignment operator.
|
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.
|
std::enable_if_t< Algorithms::SegmentsReductionKernels::isSegmentReductionKernel< SegmentsReductionKernel >::value > | 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.
|
std::enable_if_t< Algorithms::SegmentsReductionKernels::isSegmentReductionKernel< SegmentsReductionKernel >::value > | 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.
|
void | sequentialForAllRows (Function &&function) const |
| This method calls sequentialForRows for all matrix rows (for constant instances).
|
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.
|
void | transposedVectorProduct (const InVector &inVector, OutVector &outVector, typename ChooseSparseMatrixComputeReal< double, int >::type matrixMultiplicator=1.0, typename ChooseSparseMatrixComputeReal< double, int >::type outVectorMultiplicator=0.0, int begin=0, int end=0) const |
| Computes product of transposed matrix and vector.
|
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.
|
__cuda_callable__ | MatrixBase ()=default |
| Basic constructor with no parameters.
|
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__ const ValuesViewType & | getValues () const |
| Returns a constant reference to a vector with the matrix elements values.
|
__cuda_callable__ MatrixBase & | operator= (const MatrixBase &)=delete |
| Copy-assignment operator.
|
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
-
Real | is 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). |
Device | is a device where the matrix is allocated. |
Index | is a type for indexing of the matrix elements. |
MatrixType | specifies 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. |
Segments | is 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. |
ComputeReal | is the same as Real mostly but for binary matrices it is set to Index type. This can be changed by the user, of course. |
RealAllocator | is allocator for the matrix elements values. |
IndexAllocator | is allocator for the matrix elements column indexes. |