Template Numerical Library version\ main:ab3532f7
|
Implementation of sparse matrix, i.e. matrix storing only non-zero elements. More...
#include <TNL/Matrices/SparseMatrix.h>
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< double, Devices::Host, int, GeneralMatrix, Algorithms::Segments::CSR< Devices::Host, int, typename Allocators::Default< Devices::Host >::template Allocator< int > >::ViewType, typename ChooseSparseMatrixComputeReal< double, int >::type > | |
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. | |
Public Types inherited from TNL::Matrices::MatrixBase< double, Devices::Host, int, GeneralMatrix, Algorithms::Segments::CSR< Devices::Host, int, typename Allocators::Default< Devices::Host >::template Allocator< int > >::ViewType::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_> | |
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 | 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< double, Devices::Host, int, GeneralMatrix, Algorithms::Segments::CSR< Devices::Host, int, typename Allocators::Default< Devices::Host >::template Allocator< int > >::ViewType, typename ChooseSparseMatrixComputeReal< double, int >::type > | |
__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. | |
void | forAllElements (Function &&function) |
This method calls forElements for all matrix rows. | |
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 | forAllRows (Function &&function) const |
Method for parallel iteration over all matrix rows for constant instances. | |
void | forElements (IndexType begin, IndexType end, Function &&function) |
Method for iteration over all matrix rows for non-constant instances. | |
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) . | |
void | forRows (IndexType begin, IndexType end, Function &&function) const |
Method for parallel iteration over matrix rows from interval [begin, end) for constant instances. | |
ColumnIndexesViewType & | getColumnIndexes () |
Getter of column indexes for nonconstant instances. | |
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 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. | |
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. | |
const SegmentsViewType & | getSegments () const |
Getter of segments for constant instances. | |
bool | operator!= (const Matrix &matrix) const |
Comparison operator with another arbitrary matrix type. | |
SparseMatrixBase & | operator= (const SparseMatrixBase &)=delete |
Copy-assignment operator. | |
SparseMatrixBase & | operator= (SparseMatrixBase &&)=delete |
Move-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 > | reduceAllRows (Fetch &&fetch, const Reduce &reduce, Keep &&keep, const SegmentsReductionKernel &kernel=SegmentsReductionKernel{}) const |
Method for performing general reduction on all matrix rows for constant instances with function object instead of lambda function for reduction. | |
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. | |
std::enable_if_t< Algorithms::SegmentsReductionKernels::isSegmentReductionKernel< SegmentsReductionKernel >::value > | reduceRows (IndexType begin, IndexType end, Fetch &&fetch, const Reduce &reduce, Keep &&keep, const SegmentsReductionKernel &kernel=SegmentsReductionKernel{}) const |
Method for performing general reduction on matrix rows for constant instances with function object instead of lamda function for reduction. | |
void | sequentialForAllRows (Function &&function) |
This method calls sequentialForRows for all matrix rows. | |
void | sequentialForAllRows (Function &&function) const |
This method calls sequentialForRows for all matrix rows (for constant instances). | |
void | sequentialForRows (IndexType begin, IndexType end, Function &&function) |
Method for sequential iteration over all matrix rows for non-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. | |
void | vectorProduct (const InVector &inVector, OutVector &outVector, const SegmentsReductionKernel &kernel) const |
Public Member Functions inherited from TNL::Matrices::MatrixBase< double, Devices::Host, int, GeneralMatrix, Algorithms::Segments::CSR< Devices::Host, int, typename Allocators::Default< Devices::Host >::template Allocator< int > >::ViewType::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__ ValuesViewType & | getValues () |
Returns a reference to a vector with the matrix elements values. | |
__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. | |
__cuda_callable__ MatrixBase & | operator= (MatrixBase &&)=delete |
Move-assignment operator. | |
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< double, Devices::Host, int, GeneralMatrix, Algorithms::Segments::CSR< Devices::Host, int, typename Allocators::Default< Devices::Host >::template Allocator< int > >::ViewType, typename ChooseSparseMatrixComputeReal< double, int >::type > | |
ColumnIndexesViewType | columnIndexes |
SegmentsViewType | segments |
Protected Attributes inherited from TNL::Matrices::MatrixBase< double, Devices::Host, int, GeneralMatrix, Algorithms::Segments::CSR< Devices::Host, int, typename Allocators::Default< Devices::Host >::template Allocator< int > >::ViewType::getOrganization() > | |
IndexType | columns |
IndexType | rows |
ValuesViewType | values |
Additional Inherited Members | |
Protected Member Functions inherited from TNL::Matrices::SparseMatrixBase< double, Devices::Host, int, GeneralMatrix, Algorithms::Segments::CSR< Devices::Host, int, typename Allocators::Default< Devices::Host >::template Allocator< int > >::ViewType, typename ChooseSparseMatrixComputeReal< double, int >::type > | |
__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< double, Devices::Host, int, GeneralMatrix, Algorithms::Segments::CSR< Devices::Host, int, typename Allocators::Default< Devices::Host >::template Allocator< int > >::ViewType::getOrganization() > | |
__cuda_callable__ void | bind (IndexType rows, IndexType columns, ValuesViewType values) |
Re-initializes the internal attributes of the base class. | |
Implementation of sparse matrix, i.e. matrix storing only non-zero elements.
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. |
using TNL::Matrices::SparseMatrix< Real, Device, Index, MatrixType_, Segments, ComputeReal, RealAllocator, IndexAllocator >::ConstViewType |
Matrix view type for constant instances.
See SparseMatrixView.
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.
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.
realAllocator | is used for allocation of matrix elements values. |
indexAllocator | is used for allocation of matrix elements column indexes. |
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.
rows | is number of matrix rows. |
columns | is number of matrix columns. |
realAllocator | is used for allocation of matrix elements values. |
indexAllocator | is used for allocation of matrix elements column indexes. |
|
explicit |
Constructor with matrix rows capacities and number of columns.
The number of matrix rows is given by the size of rowCapacities list.
ListIndex | is the initializer list values type. |
rowCapacities | is a list telling how many matrix elements must be allocated in each row. |
columns | is the number of matrix columns. |
realAllocator | is used for allocation of matrix elements values. |
indexAllocator | is used for allocation of matrix elements column indexes. |
|
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.
RowCapacitiesVector | is the row capacities vector type. Usually it is some of TNL::Containers::Array, TNL::Containers::ArrayView, TNL::Containers::Vector or TNL::Containers::VectorView. |
rowCapacities | is a vector telling how many matrix elements must be allocated in each row. |
columns | is the number of matrix columns. |
realAllocator | is used for allocation of matrix elements values. |
indexAllocator | is used for allocation of matrix elements column indexes. |
|
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 }, ... }.
rows | is number of matrix rows. |
columns | is number of matrix columns. |
data | is a list of matrix elements values. |
realAllocator | is used for allocation of matrix elements values. |
encoding | defines encoding for sparse symmetric matrices - see TNL::Matrices::SymmetricMatrixEncoding. |
indexAllocator | is used for allocation of matrix elements column indexes. |
|
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.
MapIndex | is a type for indexing rows and columns. |
MapValue | is a type for matrix elements values in the map. |
rows | is number of matrix rows. |
columns | is number of matrix columns. |
map | is std::map containing matrix elements. |
encoding | defines encoding for sparse symmetric matrices - see TNL::Matrices::SymmetricMatrixEncoding. |
realAllocator | is used for allocation of matrix elements values. |
indexAllocator | is used for allocation of matrix elements column indexes. |
|
nodiscard |
Returns a non-modifiable view of the sparse matrix.
See SparseMatrixView.
|
staticnodiscard |
Returns string with serialization type.
The string has a form Matrices::SparseMatrix< RealType, [any_device], IndexType, General/Symmetric, Format, [any_allocator] >
.
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.
Real2 | is the real type of the input matrix. |
Index2 | is the index type of the input matrix. |
Segments2 | is the type of the segments of the input matrix. |
matrix | is the input matrix. |
matrixMultiplicator | is the factor by which the matrix is multiplied. |
|
nodiscard |
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.
fileName | is name of the file. |
|
overridevirtual |
Method for loading the matrix from a file.
file | is the input file. |
Reimplemented from TNL::Object.
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.
matrix | is input matrix for the assignment. |
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.
matrix | is input matrix for the assignment. |
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.
matrix | is input matrix for the assignment. |
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.
matrix | is input matrix for the assignment. |
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.
fileName | is name of the file. |
|
overridevirtual |
Method for saving the matrix to a file.
file | is the output file. |
Reimplemented from TNL::Object.
|
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.
columns | is the number of matrix columns. |
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.
rows | is the number of matrix rows. |
columns | is the number of matrix columns. |
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 }, ... }.
data | is a initializer list of initializer lists representing list of matrix rows. |
encoding | defines encoding for sparse symmetric matrices - see TNL::Matrices::SymmetricMatrixEncoding. |
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.
MapIndex | is a type for indexing rows and columns. |
MapValue | is a type for matrix elements values in the map. |
map | is std::map containing matrix elements. |
encoding | defines encoding for sparse symmetric matrices - see TNL::Matrices::SymmetricMatrixEncoding. |
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.
Matrix | is matrix type. This can be any matrix having methods getRows and getColumns. |
matrix | in the input matrix dimensions of which are to be adopted. |
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.
RowCapacitiesVector | is a type of vector/array used for row capacities setting. |
rowCapacities | is a vector telling the number of required non-zero matrix elements in each row. |