Template Numerical Library version\ main:bb09b17
|
Implementation of sparse matrix view. More...
#include <TNL/Matrices/SparseMatrixBase.h>
Public Types | |
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 | |
__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. | |
ColumnIndexesViewType & | getColumnIndexes () |
Getter of column indexes for nonconstant instances. | |
const ColumnIndexesViewType & | getColumnIndexes () 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. | |
SegmentsViewType & | getSegments () |
Getter of segments for non-constant instances. | |
const SegmentsViewType & | getSegments () const |
Getter of segments for constant instances. | |
template<typename Matrix > | |
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. | |
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__ 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. | |
bool | operator!= (const Matrix &matrix) const |
Comparison operator with another arbitrary matrix view type. | |
bool | operator!= (const MatrixT &matrix) const |
__cuda_callable__ MatrixBase & | operator= (const MatrixBase &)=delete |
Copy-assignment operator. | |
__cuda_callable__ MatrixBase & | operator= (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::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 | isSymmetric () |
Test of symmetric matrix type. | |
Protected Member Functions | |
__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. | |
Protected Attributes | |
ColumnIndexesViewType | columnIndexes |
SegmentsViewType | segments |
Protected Attributes inherited from TNL::Matrices::MatrixBase< Real, Device, Index, MatrixType, SegmentsView::getOrganization() > | |
IndexType | columns |
IndexType | rows |
ValuesViewType | values |
Implementation of sparse matrix view.
It serves as an accessor to SparseMatrix for example when passing the matrix to lambda functions. SparseMatrix view can be also created in CUDA kernels.
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. |
using TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::ColumnIndexesViewType |
using TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::DefaultSegmentsReductionKernel = typename Algorithms::SegmentsReductionKernels::DefaultKernel< SegmentsView >::type |
Type of the kernel used for parallel reductions on segments.
We are assuming that the default segments reduction kernel provides a static reduceAllSegments method, and thus it does not have to be instantiated and initialized. If the user wants to use a more complicated kernel, such as CSRAdaptive, it must be instantiated and initialized by the user and the object must be passed to the vectorProduct or reduceRows method.
using TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::RowView |
Type for accessing matrix rows.
__cuda_callable__ TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::SparseMatrixBase | ( | IndexType | rows, |
IndexType | columns, | ||
typename Base::ValuesViewType | values, | ||
ColumnIndexesViewType | columnIndexes, | ||
SegmentsViewType | segments ) |
Constructor with all necessary data and views.
rows | is a number of matrix rows. |
columns | is a number of matrix columns. |
values | is a vector view with matrix elements values. |
columnIndexes | is a vector view with matrix elements column indexes. |
segments | is a segments view representing the sparse matrix format. |
|
default |
Copy constructor.
matrix | is an input sparse matrix view. |
|
defaultnoexcept |
Move constructor.
matrix | is an input sparse matrix view. |
__cuda_callable__ void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::addElement | ( | IndexType | row, |
IndexType | column, | ||
const RealType & | value, | ||
const RealType & | thisElementMultiplicator = 1.0 ) |
Add element at given row and column to given value.
This method can be called from the host system (CPU) no matter where the matrix is allocated. If the matrix is allocated on GPU this method can be called even from device kernels. If the matrix is allocated in GPU device this method is called from CPU, it transfers values of each matrix element separately and so the performance is very low. For higher performance see. SparseMatrix::getRow or SparseMatrix::forElements and SparseMatrix::forAllElements. The call may fail if the matrix row capacity is exhausted.
row | is row index of the element. |
column | is columns index of the element. |
value | is the value the element will be set to. |
thisElementMultiplicator | is multiplicator the original matrix element value is multiplied by before addition of given value. |
|
protected |
Re-initializes the internal attributes of the base class.
Note that this function is protected to ensure that the user cannot modify the base class of a matrix. For the same reason, in future code development we also need to make sure that all non-const functions in the base class return by value and not by reference.
|
nodiscard |
Finds element in the matrix and returns its position in the arrays with values and columnIndexes.
If the element is not found, the method returns the padding index.
row | is the row index of the element. |
column | is the column index of the element. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::forAllElements | ( | Function && | function | ) |
This method calls forElements for all matrix rows.
See SparseMatrix::forElements.
Function | is a type of lambda function that will operate on matrix elements. |
function | is an instance of the lambda function to be called in each row. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::forAllElements | ( | Function && | function | ) | const |
This method calls forElements for all matrix rows (for constant instances).
See SparseMatrix::forElements.
Function | is a type of lambda function that will operate on matrix elements. |
function | is an instance of the lambda function to be called in each row. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::forAllRows | ( | Function && | function | ) |
Method for parallel iteration over all matrix rows.
In each row, given lambda function is performed. Each row is processed by at most one thread unlike the method SparseMatrixBase::forAllElements where more than one thread can be mapped to each row.
Function | is type of the lambda function. |
function | is an instance of the lambda function to be called for each row. |
RowView represents matrix row - see TNL::Matrices::SparseMatrixBase::RowView.
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::forAllRows | ( | Function && | function | ) | const |
Method for parallel iteration over all matrix rows for constant instances.
In each row, given lambda function is performed. Each row is processed by at most one thread unlike the method SparseMatrixBase::forAllElements where more than one thread can be mapped to each row.
Function | is type of the lambda function. |
function | is an instance of the lambda function to be called for each row. |
ConstRowView represents matrix row - see TNL::Matrices::SparseMatrixBase::ConstRowView.
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::forElements | ( | IndexType | begin, |
IndexType | end, | ||
Function && | function ) |
Method for iteration over all matrix rows for non-constant instances.
Function | is type of lambda function that will operate on matrix elements. It should have form like |
The localIdx parameter is a rank of the non-zero element in given row.
begin | defines beginning of the range [begin, end) of rows to be processed. |
end | defines ending of the range [begin, end) of rows to be processed. |
function | is an instance of the lambda function to be called in each row. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::forElements | ( | IndexType | begin, |
IndexType | end, | ||
Function && | function ) const |
Method for iteration over all matrix rows for constant instances.
Function | is type of lambda function that will operate on matrix elements. It should have form like |
The localIdx parameter is a rank of the non-zero element in given row.
begin | defines beginning of the range [begin, end) of rows to be processed. |
end | defines ending of the range [begin, end) of rows to be processed. |
function | is an instance of the lambda function to be called in each row. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::forRows | ( | IndexType | begin, |
IndexType | end, | ||
Function && | function ) |
Method for parallel iteration over matrix rows from interval [begin, end)
.
In each row, given lambda function is performed. Each row is processed by at most one thread unlike the method SparseMatrixBase::forElements where more than one thread can be mapped to each row.
Function | is type of the lambda function. |
begin | defines beginning of the range [begin, end) of rows to be processed. |
end | defines ending of the range [begin, end) of rows to be processed. |
function | is an instance of the lambda function to be called for each row. |
RowView represents matrix row - see TNL::Matrices::SparseMatrixBase::RowView.
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::forRows | ( | IndexType | begin, |
IndexType | end, | ||
Function && | function ) const |
Method for parallel iteration over matrix rows from interval [begin, end)
for constant instances.
In each row, given lambda function is performed. Each row is processed by at most one thread unlike the method SparseMatrixBase::forElements where more than one thread can be mapped to each row.
Function | is type of the lambda function. |
begin | defines beginning of the range [begin, end) of rows to be processed. |
end | defines ending of the range [begin, end) of rows to be processed. |
function | is an instance of the lambda function to be called for each row. |
ConstRowView represents matrix row - see TNL::Matrices::SparseMatrixBase::ConstRowView.
|
nodiscard |
Getter of column indexes for nonconstant instances.
|
nodiscard |
Getter of column indexes for constant instances.
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::getCompressedRowLengths | ( | Vector & | rowLengths | ) | const |
Computes number of non-zeros in each row.
rowLengths | is a vector into which the number of non-zeros in each row will be stored. |
|
nodiscard |
Returns value of matrix element at position given by its row and column index.
This method can be called from the host system (CPU) no matter where the matrix is allocated. If the matrix is allocated on GPU this method can be called even from device kernels. If the matrix is allocated in GPU device this method is called from CPU, it transfers values of each matrix element separately and so the performance is very low. For higher performance see. SparseMatrix::getRow or SparseMatrix::forElements and SparseMatrix::forAllElements.
row | is a row index of the matrix element. |
column | i a column index of the matrix element. |
|
nodiscardoverridevirtual |
Returns number of non-zero matrix elements.
This method really counts the non-zero matrix elements and so it returns zero for matrix having all allocated elements set to zero.
Reimplemented from TNL::Matrices::MatrixBase< Real, Device, Index, MatrixType, SegmentsView::getOrganization() >.
|
nodiscard |
Non-constant getter of simple structure for accessing given matrix row.
rowIdx | is matrix row index. |
See SparseMatrixRowView.
|
nodiscard |
Constant getter of simple structure for accessing given matrix row.
rowIdx | is matrix row index. |
See SparseMatrixRowView.
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::getRowCapacities | ( | Vector & | rowCapacities | ) | const |
Compute capacities of all rows.
The row capacities are not stored explicitly and must be computed.
rowCapacities | is a vector where the row capacities will be stored. |
|
nodiscard |
Returns capacity of given matrix row.
row | index of matrix row. |
|
nodiscard |
Getter of segments for non-constant instances.
Segments are a structure for addressing the matrix elements columns and values. In fact, Segments represent the sparse matrix format.
|
nodiscard |
Getter of segments for constant instances.
Segments are a structure for addressing the matrix elements columns and values. In fact, Segments represent the sparse matrix format.
|
staticnodiscard |
Returns string with serialization type.
The string has a form Matrices::SparseMatrix< RealType, [any_device], IndexType, General/Symmetric, Format, [any_allocator] >
.
|
nodiscard |
Comparison operator with another arbitrary matrix type.
matrix | is the right-hand side matrix. |
|
delete |
Copy-assignment operator.
It is a deleted function, because sparse matrix assignment in general requires reallocation.
|
nodiscard |
Comparison operator with another arbitrary matrix type.
matrix | is the right-hand side matrix. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::print | ( | std::ostream & | str | ) | const |
Method for printing the matrix to output stream.
str | is the output stream. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::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.
Fetch | is a type of lambda function for data fetch declared as |
The return type of this lambda can be any non void.
Reduce | is a type of lambda function for reduction declared as |
Keep | is a type of lambda function for storing results of reduction in each row. It is declared as |
FetchValue | is type returned by the Fetch lambda function. |
fetch | is an instance of lambda function for data fetch. |
reduce | is an instance of lambda function for reduction. |
keep | in an instance of lambda function for storing results. |
identity | is the identity element for the reduction operation, i.e. element which does not change the result of the reduction. |
kernel | is an instance of the segments reduction kernel to be used for the operation. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::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.
Fetch | is a type of lambda function for data fetch declared as |
The return type of this lambda can be any non void.
Reduce | is a type of lambda function for reduction declared as |
Keep | is a type of lambda function for storing results of reduction in each row. It is declared as |
FetchValue | is type returned by the Fetch lambda function. |
begin | defines beginning of the range [begin, end) of rows to be processed. |
end | defines ending of the range [begin, end) of rows to be processed. |
fetch | is an instance of lambda function for data fetch. |
reduce | is an instance of lambda function for reduction. |
keep | in an instance of lambda function for storing results. |
identity | is the identity element for the reduction operation, i.e. element which does not change the result of the reduction. |
kernel | is an instance of the segments reduction kernel to be used for the operation. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::sequentialForAllRows | ( | Function && | function | ) |
This method calls sequentialForRows for all matrix rows.
See SparseMatrixBase::sequentialForAllRows.
Function | is a type of lambda function that will operate on matrix elements. |
function | is an instance of the lambda function to be called in each row. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::sequentialForAllRows | ( | Function && | function | ) | const |
This method calls sequentialForRows for all matrix rows (for constant instances).
See SparseMatrixBase::sequentialForRows.
Function | is a type of lambda function that will operate on matrix elements. |
function | is an instance of the lambda function to be called in each row. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::sequentialForRows | ( | IndexType | begin, |
IndexType | end, | ||
Function && | function ) |
Method for sequential iteration over all matrix rows for non-constant instances.
Function | is type of lambda function that will operate on matrix elements. It should have form like |
RowView represents matrix row - see TNL::Matrices::SparseMatrixBase::RowView.
begin | defines beginning of the range [begin, end) of rows to be processed. |
end | defines ending of the range [begin, end) of rows to be processed. |
function | is an instance of the lambda function to be called in each row. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::sequentialForRows | ( | IndexType | begin, |
IndexType | end, | ||
Function && | function ) const |
Method for sequential iteration over all matrix rows for constant instances.
Function | is type of lambda function that will operate on matrix elements. It should have form like |
ConstRowView represents matrix row - see TNL::Matrices::SparseMatrixBase::ConstRowView.
begin | defines beginning of the range [begin, end) of rows to be processed. |
end | defines ending of the range [begin, end) of rows to be processed. |
function | is an instance of the lambda function to be called in each row. |
__cuda_callable__ void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::setElement | ( | IndexType | row, |
IndexType | column, | ||
const RealType & | value ) |
Sets element at given row and column to given value.
This method can be called from the host system (CPU) no matter where the matrix is allocated. If the matrix is allocated on GPU this method can be called even from device kernels. If the matrix is allocated in GPU device this method is called from CPU, it transfers values of each matrix element separately and so the performance is very low. For higher performance see. SparseMatrix::getRow or SparseMatrix::forElements and SparseMatrix::forAllElements. The call may fail if the matrix row capacity is exhausted.
row | is row index of the element. |
column | is columns index of the element. |
value | is the value the element will be set to. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::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.
InVector | is type of input vector. It can be TNL::Containers::Vector, TNL::Containers::VectorView, TNL::Containers::Array, TNL::Containers::ArrayView, or similar container. |
OutVector | is type of output vector. It can be TNL::Containers::Vector, TNL::Containers::VectorView, TNL::Containers::Array, TNL::Containers::ArrayView, or similar container. |
inVector | is input vector. |
outVector | is output vector. |
matrixMultiplicator | is a factor by which the matrix is multiplied. It is one by default. |
outVectorMultiplicator | is a factor by which the outVector is multiplied before added |
begin | is the beginning of the rows range for which the vector product is computed. It is zero by default. |
end | is the end of the rows range for which the vector product is computed. It is number if the matrix rows by default. |
void TNL::Matrices::SparseMatrixBase< Real, Device, Index, MatrixType, SegmentsView, ComputeReal >::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.
More precisely, it computes:
InVector | is type of input vector. It can be TNL::Containers::Vector, TNL::Containers::VectorView, TNL::Containers::Array, TNL::Containers::ArrayView, or similar container. |
OutVector | is type of output vector. It can be TNL::Containers::Vector, TNL::Containers::VectorView, TNL::Containers::Array, TNL::Containers::ArrayView, or similar container. |
inVector | is input vector. |
outVector | is output vector. |
matrixMultiplicator | is a factor by which the matrix is multiplied. It is one by default. |
outVectorMultiplicator | is a factor by which the outVector is multiplied before added to the result of matrix-vector product. It is zero by default. |
begin | is the beginning of the rows range for which the vector product is computed. It is zero by default. |
end | is the end of the rows range for which the vector product is computed. It is number if the matrix rows by default. |
kernel | is an instance of the segments reduction kernel to be used for the operation. |