Template Numerical Library version\ main:94209208
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization > Class Template Reference

A common base class for TridiagonalMatrix and TridiagonalMatrixView. More...

#include <TNL/Matrices/TridiagonalMatrixBase.h>

Inheritance diagram for TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >:
Inheritance graph
[legend]
Collaboration diagram for TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >:
Collaboration graph
[legend]

Public Types

using ConstRowView = typename RowView::ConstRowView
 Type for accessing constant matrix rows.
 
using DeviceType = Device
 The device where the matrix is allocated.
 
using IndexerType = details::TridiagonalMatrixIndexer< Index, Organization >
 
using IndexType = Index
 The type used for matrix elements indexing.
 
using RealType = typename Base::RealType
 The type of matrix elements.
 
using RowView = TridiagonalMatrixRowView< typename Base::ValuesViewType, IndexerType >
 Type for accessing matrix rows.
 
- Public Types inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
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__ TridiagonalMatrixBase ()=default
 Constructor with no parameters.
 
__cuda_callable__ TridiagonalMatrixBase (const TridiagonalMatrixBase &)=default
 Copy constructor.
 
__cuda_callable__ TridiagonalMatrixBase (TridiagonalMatrixBase &&) noexcept=default
 Move constructor.
 
__cuda_callable__ TridiagonalMatrixBase (typename Base::ValuesViewType values, IndexerType indexer)
 Constructor with all necessary data and views.
 
__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.
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_>
void addMatrix (const TridiagonalMatrixBase< Real_, Device_, Index_, Organization_ > &matrix, const RealType &matrixMultiplicator=1.0, const RealType &thisMatrixMultiplicator=1.0)
 
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.
 
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.
 
__cuda_callable__ IndexerTypegetIndexer ()
 This method returns matrix elements indexer used by this matrix.
 
__cuda_callable__ const IndexerTypegetIndexer () const
 This method returns matrix elements indexer used by this matrix.
 
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.
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_>
bool operator!= (const TridiagonalMatrixBase< Real_, Device_, Index_, Organization_ > &matrix) const
 Comparison operator with another multidiagonal matrix.
 
TridiagonalMatrixBaseoperator= (const TridiagonalMatrixBase &)=delete
 Copy-assignment operator.
 
TridiagonalMatrixBaseoperator= (TridiagonalMatrixBase &&)=delete
 Move-assignment operator.
 
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_>
bool operator== (const TridiagonalMatrixBase< Real_, Device_, Index_, Organization_ > &matrix) const
 Comparison operator with another tridiagonal matrix.
 
void print (std::ostream &str) const
 Method for printing the matrix to output stream.
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void reduceAllRows (Fetch &fetch, Reduce &reduce, Keep &keep, const FetchReal &identity) const
 Method for performing general reduction on all matrix rows for constant instances.
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void reduceRows (IndexType begin, IndexType end, Fetch &fetch, Reduce &reduce, Keep &keep, const FetchReal &identity) 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 setValue (const RealType &value)
 Set all matrix elements to given value.
 
template<typename InVector , typename OutVector >
void vectorProduct (const InVector &inVector, OutVector &outVector, RealType matrixMultiplicator=1.0, RealType outVectorMultiplicator=0.0, IndexType begin=0, IndexType end=0) const
 Computes product of matrix and vector.
 
- Public Member Functions inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
__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::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
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 Member Functions

__cuda_callable__ void bind (typename Base::ValuesViewType values, IndexerType indexer)
 Re-initializes the internal attributes of the base class.
 
__cuda_callable__ IndexType getElementIndex (IndexType row, IndexType column) const
 
- Protected Member Functions inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
__cuda_callable__ void bind (IndexType rows, IndexType columns, ValuesViewType values)
 Re-initializes the internal attributes of the base class.
 

Protected Attributes

IndexerType indexer
 
- Protected Attributes inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
IndexType columns
 
IndexType rows
 
ValuesViewType values
 

Detailed Description

template<typename Real, typename Device, typename Index, ElementsOrganization Organization>
class TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >

A common base class for TridiagonalMatrix and TridiagonalMatrixView.

Template Parameters
Realis a type of matrix elements.
Deviceis a device where the matrix is allocated.
Indexis a type for indexing of the matrix elements.
Organizationtells the ordering of matrix elements. It is either RowMajorOrder or ColumnMajorOrder.

Constructor & Destructor Documentation

◆ TridiagonalMatrixBase()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
__cuda_callable__ TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::TridiagonalMatrixBase ( typename Base::ValuesViewType values,
IndexerType indexer )

Constructor with all necessary data and views.

Parameters
valuesis a vector view with matrix elements values
indexeris an indexer of matrix elements

Member Function Documentation

◆ addElement()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
__cuda_callable__ void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::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. TridiagonalMatrix::getRow or TridiagonalMatrix::forElements and TridiagonalMatrix::forAllElements. The call may fail if the matrix row capacity is exhausted.

Parameters
rowis row index of the element.
columnis columns index of the element.
valueis the value the element will be set to.
thisElementMultiplicatoris multiplicator the original matrix element value is multiplied by before addition of given value.
Example
#include <iostream>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void
addElements()
{
const int matrixSize( 5 );
TNL::Matrices::TridiagonalMatrix< double, Device > matrix( matrixSize, // number of rows
matrixSize // number of columns
);
auto view = matrix.getView();
for( int i = 0; i < matrixSize; i++ )
view.setElement( i, i, i ); // or matrix.setElement
std::cout << "Initial matrix is: " << std::endl << matrix << std::endl;
for( int i = 0; i < matrixSize; i++ ) {
if( i > 0 )
view.addElement( i, i - 1, 1.0, 5.0 ); // or matrix.addElement
view.addElement( i, i, 1.0, 5.0 ); // or matrix.addElement
if( i < matrixSize - 1 )
view.addElement( i, i + 1, 1.0, 5.0 ); // or matrix.addElement
}
std::cout << "Matrix after addition is: " << std::endl << matrix << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Add elements on host:" << std::endl;
addElements< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Add elements on CUDA device:" << std::endl;
addElements< TNL::Devices::Cuda >();
#endif
}
Implementation of sparse tridiagonal matrix.
Definition TridiagonalMatrix.h:56
T endl(T... args)
Output
Add elements on host:
Initial matrix is:
Row: 0 ->
Row: 1 -> 1:1
Row: 2 -> 2:2
Row: 3 -> 3:3
Row: 4 -> 4:4
Matrix after addition is:
Row: 0 -> 0:1 1:1
Row: 1 -> 0:1 1:6 2:1
Row: 2 -> 1:1 2:11 3:1
Row: 3 -> 2:1 3:16 4:1
Row: 4 -> 3:1 4:21
Add elements on CUDA device:
Initial matrix is:
Row: 0 ->
Row: 1 -> 1:1
Row: 2 -> 2:2
Row: 3 -> 3:3
Row: 4 -> 4:4
Matrix after addition is:
Row: 0 -> 0:1 1:1
Row: 1 -> 0:1 1:6 2:1
Row: 2 -> 1:1 2:11 3:1
Row: 3 -> 2:1 3:16 4:1
Row: 4 -> 3:1 4:21

◆ bind()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
__cuda_callable__ void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::bind ( typename Base::ValuesViewType values,
IndexerType indexer )
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.

◆ forAllElements() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Function >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::forAllElements ( Function & function)

This method calls forElements for all matrix rows.

See TridiagonalMatrix::forElements.

Template Parameters
Functionis a type of lambda function that will operate on matrix elements.
Parameters
functionis an instance of the lambda function to be called in each row.
Example
#include <iostream>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
forAllElementsExample()
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 1 3 . . . \ -> { 0, 1, 3 }
* | 2 1 3 . . | -> { 2, 1, 3 }
* | . 2 1 3 . | -> { 2, 1, 3 }
* | . . 2 1 3 | -> { 2, 1, 3 }
* \ . . . 2 1 / 0 -> { 2, 1, 0 }
*/
TNL::Matrices::TridiagonalMatrix< double, Device > matrix( 5, // number of matrix rows
5 ); // number of matrix columns
auto view = matrix.getView();
auto f = [] __cuda_callable__( int rowIdx, int localIdx, int columnIdx, double& value )
{
/***
* 'forElements' method iterates only over matrix elements lying on given subdiagonals
* and so we do not need to check anything. The element value can be expressed
* by the 'localIdx' variable, see the following figure:
*
* 0 1 2 <- localIdx values
* -------
* 0 / 1 3 . . . \ -> { 0, 1, 3 }
* | 2 1 3 . . | -> { 2, 1, 3 }
* | . 2 1 3 . | -> { 2, 1, 3 }
* | . . 2 1 3 | -> { 2, 1, 3 }
* \ . . . 2 1 / 0 -> { 2, 1, 0 }
*
*/
value = 3 - localIdx;
};
view.forAllElements( f ); // or matrix.forAllElements
std::cout << matrix << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forAllElementsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forAllElementsExample< TNL::Devices::Cuda >();
#endif
}
#define __cuda_callable__
Definition Macros.h:49
Output
Creating matrix on host:
Row: 0 -> 0:2 1:1
Row: 1 -> 0:3 1:2 2:1
Row: 2 -> 1:3 2:2 3:1
Row: 3 -> 2:3 3:2 4:1
Row: 4 -> 3:3 4:2
Creating matrix on CUDA device:
Row: 0 -> 0:2 1:1
Row: 1 -> 0:3 1:2 2:1
Row: 2 -> 1:3 2:2 3:1
Row: 3 -> 2:3 3:2 4:1
Row: 4 -> 3:3 4:2

◆ forAllElements() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Function >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::forAllElements ( Function & function) const

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

See TridiagonalMatrix::forElements.

Template Parameters
Functionis a type of lambda function that will operate on matrix elements.
Parameters
functionis an instance of the lambda function to be called in each row.
Example
#include <iostream>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
forAllElementsExample()
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 1 3 . . . \ -> { 0, 1, 3 }
* | 2 1 3 . . | -> { 2, 1, 3 }
* | . 2 1 3 . | -> { 2, 1, 3 }
* | . . 2 1 3 | -> { 2, 1, 3 }
* \ . . . 2 1 / 0 -> { 2, 1, 0 }
*/
TNL::Matrices::TridiagonalMatrix< double, Device > matrix( 5, // number of matrix rows
5 ); // number of matrix columns
auto view = matrix.getView();
auto f = [] __cuda_callable__( int rowIdx, int localIdx, int columnIdx, double& value )
{
/***
* 'forElements' method iterates only over matrix elements lying on given subdiagonals
* and so we do not need to check anything. The element value can be expressed
* by the 'localIdx' variable, see the following figure:
*
* 0 1 2 <- localIdx values
* -------
* 0 / 1 3 . . . \ -> { 0, 1, 3 }
* | 2 1 3 . . | -> { 2, 1, 3 }
* | . 2 1 3 . | -> { 2, 1, 3 }
* | . . 2 1 3 | -> { 2, 1, 3 }
* \ . . . 2 1 / 0 -> { 2, 1, 0 }
*
*/
value = 3 - localIdx;
};
view.forAllElements( f ); // or matrix.forAllElements
std::cout << matrix << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forAllElementsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forAllElementsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:2 1:1
Row: 1 -> 0:3 1:2 2:1
Row: 2 -> 1:3 2:2 3:1
Row: 3 -> 2:3 3:2 4:1
Row: 4 -> 3:3 4:2
Creating matrix on CUDA device:
Row: 0 -> 0:2 1:1
Row: 1 -> 0:3 1:2 2:1
Row: 2 -> 1:3 2:2 3:1
Row: 3 -> 2:3 3:2 4:1
Row: 4 -> 3:3 4:2

◆ forAllRows() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Function >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::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 TridiagonalMatrixBase::forAllElements where more than one thread can be mapped to each row.

Template Parameters
Functionis type of the lambda function.
Parameters
functionis an instance of the lambda function to be called for each row.
auto function = [] __cuda_callable__ ( RowView& row ) { ... };
RowView is a simple structure for accessing rows of tridiagonal matrix.
Definition TridiagonalMatrixRowView.h:27

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

Example
#include <iostream>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
forRowsExample()
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 2 . . . . \ -> { 0, 0, 1 }
* | -1 2 -1 . . | -> { 0, 2, 1 }
* | . -1 2 -1. . | -> { 3, 2, 1 }
* | . . -1 2 -1 | -> { 3, 2, 1 }
* \ . . . . 2 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -1, 0, 1 }.
*/
const int size = 5;
MatrixType matrix( size, size );
auto view = matrix.getView();
auto f = [] __cuda_callable__( typename MatrixType::RowView & row )
{
const int& rowIdx = row.getRowIndex();
if( rowIdx > 0 )
row.setElement( 0, -1.0 ); // elements below the diagonal
row.setElement( 1, 2.0 ); // elements on the diagonal
if( rowIdx < size - 1 ) // elements above the diagonal
row.setElement( 2, -1.0 );
};
view.forAllRows( f ); // or matrix.forAllRows
std::cout << matrix << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forRowsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forRowsExample< TNL::Devices::Cuda >();
#endif
}
Structure for specifying type of sparse matrix.
Definition MatrixType.h:17
Output
Creating matrix on host:
Row: 0 -> 0:2 1:-1
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 -> 3:-1 4:2
Creating matrix on CUDA device:
Row: 0 -> 0:2 1:-1
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 -> 3:-1 4:2

◆ forAllRows() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Function >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::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 TridiagonalMatrixBase::forAllElements where more than one thread can be mapped to each row.

Template Parameters
Functionis type of the lambda function.
Parameters
functionis an instance of the lambda function to be called for each row.
auto function = [] __cuda_callable__ ( const ConstRowView& row ) { ... };
typename RowView::ConstRowView ConstRowView
Type for accessing constant matrix rows.
Definition TridiagonalMatrixBase.h:53

ConstRowView represents matrix row - see TNL::Matrices::TridiagonalMatrixBase::ConstRowView.

◆ forElements() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Function >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::forElements ( IndexType begin,
IndexType end,
Function & function )

Method for iteration over all matrix rows for non-constant instances.

Template Parameters
Functionis type of lambda function that will operate on matrix elements. It should have form like
auto function = [=] __cuda_callable__ ( IndexType rowIdx, IndexType localIdx, IndexType columnIdx, RealType& value )
{ ... };
typename Base::RealType RealType
The type of matrix elements.
Definition TridiagonalMatrixBase.h:30
Index IndexType
The type used for matrix elements indexing.
Definition TridiagonalMatrixBase.h:40

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

Parameters
begindefines beginning of the range [begin, end) of rows to be processed.
enddefines ending of the range [begin, end) of rows to be processed.
functionis an instance of the lambda function to be called in each row.
Example
#include <iostream>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
forElementsExample()
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 2 1 . . . \ -> { 0, 2, 1 }
* | 3 2 1 . . | -> { 3, 2, 1 }
* | . 3 2 1 . | -> { 3, 2, 1 }
* | . . 3 2 1 | -> { 3, 2, 1 }
* \ . . . 3 2 / 0 -> { 3, 2, 0 }
*/
TNL::Matrices::TridiagonalMatrix< double, Device > matrix( 5, // number of matrix rows
5 ); // number of matrix columns
auto view = matrix.getView();
auto f = [] __cuda_callable__( int rowIdx, int localIdx, int columnIdx, double& value )
{
/***
* 'forElements' method iterates only over matrix elements lying on given subdiagonals
* and so we do not need to check anything. The element value can be expressed
* by the 'localIdx' variable, see the following figure:
*
* 0 1 2 <- localIdx values
* -------
* 0 / 2 1 . . . \ -> { 0, 2, 1 }
* | 3 2 1 . . | -> { 3, 2, 1 }
* | . 3 2 1 . | -> { 3, 2, 1 }
* | . . 3 2 1 | -> { 3, 2, 1 }
* \ . . . 3 2 / 0 -> { 3, 2, 0 }
*
*/
value = 3 - localIdx;
};
view.forElements( 0, matrix.getRows(), f ); // or matrix.forElements
std::cout << matrix << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forElementsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forElementsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:2 1:1
Row: 1 -> 0:3 1:2 2:1
Row: 2 -> 1:3 2:2 3:1
Row: 3 -> 2:3 3:2 4:1
Row: 4 -> 3:3 4:2
Creating matrix on CUDA device:
Row: 0 -> 0:2 1:1
Row: 1 -> 0:3 1:2 2:1
Row: 2 -> 1:3 2:2 3:1
Row: 3 -> 2:3 3:2 4:1
Row: 4 -> 3:3 4:2

◆ forElements() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Function >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::forElements ( IndexType begin,
IndexType end,
Function & function ) const

Method for iteration over all matrix rows for constant instances.

Template Parameters
Functionis type of lambda function that will operate on matrix elements. It should have form like
auto function = [] __cuda_callable__ ( IndexType rowIdx, IndexType localIdx, IndexType columnIdx, const RealType& value )
{ ... };

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

Parameters
begindefines beginning of the range [begin, end) of rows to be processed.
enddefines ending of the range [begin, end) of rows to be processed.
functionis an instance of the lambda function to be called in each row.

◆ forRows() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Function >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::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 TridiagonalMatrixBase::forElements where more than one thread can be mapped to each row.

Template Parameters
Functionis type of the lambda function.
Parameters
begindefines beginning of the range [begin, end) of rows to be processed.
enddefines ending of the range [begin, end) of rows to be processed.
functionis an instance of the lambda function to be called for each row.
auto function = [] __cuda_callable__ ( RowView& row ) { ... };

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

Example
#include <iostream>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
forRowsExample()
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 2 . . . . \ -> { 0, 0, 1 }
* | -1 2 -1 . . | -> { 0, 2, 1 }
* | . -1 2 -1. . | -> { 3, 2, 1 }
* | . . -1 2 -1 | -> { 3, 2, 1 }
* \ . . . . 2 / -> { 3, 2, 1 }
*
* The diagonals offsets are { -1, 0, 1 }.
*/
const int size = 5;
MatrixType matrix( size, size );
auto view = matrix.getView();
auto f = [] __cuda_callable__( typename MatrixType::RowView & row )
{
const int& rowIdx = row.getRowIndex();
if( rowIdx > 0 )
row.setElement( 0, -1.0 ); // elements below the diagonal
row.setElement( 1, 2.0 ); // elements on the diagonal
if( rowIdx < size - 1 ) // elements above the diagonal
row.setElement( 2, -1.0 );
};
view.forAllRows( f ); // or matrix.forAllRows
std::cout << matrix << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Creating matrix on host: " << std::endl;
forRowsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrix on CUDA device: " << std::endl;
forRowsExample< TNL::Devices::Cuda >();
#endif
}
Output
Creating matrix on host:
Row: 0 -> 0:2 1:-1
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 -> 3:-1 4:2
Creating matrix on CUDA device:
Row: 0 -> 0:2 1:-1
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 -> 3:-1 4:2

◆ forRows() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Function >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::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 TridiagonalMatrixBase::forElements where more than one thread can be mapped to each row.

Template Parameters
Functionis type of the lambda function.
Parameters
begindefines beginning of the range [begin, end) of rows to be processed.
enddefines ending of the range [begin, end) of rows to be processed.
functionis an instance of the lambda function to be called for each row.
auto function = [] __cuda_callable__ ( const ConstRowView& row ) { ... };

ConstRowView represents matrix row - see TNL::Matrices::TridiagonalMatrixBase::ConstRowView.

◆ getCompressedRowLengths()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Vector >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::getCompressedRowLengths ( Vector & rowLengths) const

Computes number of non-zeros in each row.

Parameters
rowLengthsis a vector into which the number of non-zeros in each row will be stored.
Example
#include <iostream>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
laplaceOperatorMatrix()
{
const int gridSize( 6 );
const int matrixSize = gridSize;
TNL::Matrices::TridiagonalMatrix< double, Device > matrix( matrixSize, // number of rows
matrixSize // number of columns
);
matrix.setElements( {
// clang-format off
{ 0.0, 1.0 },
{ -1.0, 2.0, -1.0 },
{ -1.0, 2.0, -1.0 },
{ -1.0, 2.0, -1.0 },
{ -1.0, 2.0, -1.0 },
{ 0.0, 1.0 }
// clang-format on
} );
auto view = matrix.getView();
view.getCompressedRowLengths( rowLengths ); // or matrix.getCompressedRowLengths
std::cout << "Laplace operator matrix: " << std::endl << matrix << std::endl;
std::cout << "Compressed row lengths: " << rowLengths << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Creating Laplace operator matrix on CPU ... " << std::endl;
laplaceOperatorMatrix< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating Laplace operator matrix on CUDA GPU ... " << std::endl;
laplaceOperatorMatrix< TNL::Devices::Cuda >();
#endif
}
Vector extends Array with algebraic operations.
Definition Vector.h:36
Output
Creating Laplace operator matrix on CPU ...
Laplace operator matrix:
Row: 0 -> 0:1
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 -> 3:-1 4:2 5:-1
Row: 5 -> 5:1
Compressed row lengths: [ 1, 3, 3, 3, 3, 1 ]
Creating Laplace operator matrix on CUDA GPU ...
Laplace operator matrix:
Row: 0 -> 0:1
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 -> 3:-1 4:2 5:-1
Row: 5 -> 5:1
Compressed row lengths: [ 1, 3, 3, 3, 3, 1 ]

◆ getElement()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
__cuda_callable__ auto TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::getElement ( IndexType row,
IndexType column ) const

Returns value of matrix element at position given by its row and column index.

This method can be called from the host system (CPU) no matter where the matrix is allocated. If the matrix is allocated on GPU this method can be called even from device kernels. If the matrix is allocated in GPU device this method is called from CPU, it transfers values of each matrix element separately and so the performance is very low. For higher performance see. TridiagonalMatrix::getRow or TridiagonalMatrix::forElements and TridiagonalMatrix::forAllElements.

Parameters
rowis a row index of the matrix element.
columni a column index of the matrix element.
Returns
value of given matrix element.
Example
#include <iostream>
#include <iomanip>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void
getElements()
{
const int matrixSize( 5 );
TNL::Matrices::TridiagonalMatrix< double, Device > matrix( matrixSize, // number of matrix columns
{ // matrix elements definition
{ 0.0, 2.0, -1.0 },
{ -1.0, 2.0, -1.0 },
{ -1.0, 2.0, -1.0 },
{ -1.0, 2.0, -1.0 },
{ -1.0, 2.0, 0.0 } } );
auto view = matrix.getView();
for( int i = 0; i < matrixSize; i++ ) {
for( int j = 0; j < matrixSize; j++ )
std::cout << std::setw( 5 ) << view.getElement( i, j ); // or matrix.getElement
}
}
int
main( int argc, char* argv[] )
{
std::cout << "Get elements on host:" << std::endl;
getElements< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Get elements on CUDA device:" << std::endl;
getElements< TNL::Devices::Cuda >();
#endif
}
__cuda_callable__ RealType getElement(IndexType row, IndexType column) const
Returns value of matrix element at position given by its row and column index.
Definition TridiagonalMatrixBase.hpp:173
ViewType getView()
Returns a modifiable view of the tridiagonal matrix.
Definition TridiagonalMatrix.hpp:36
T setw(T... args)
Output
Get elements on host:
2 -1 0 0 0
-1 2 -1 0 0
0 -1 2 -1 0
0 0 -1 2 -1
0 0 0 -1 2
Get elements on CUDA device:
2 -1 0 0 0
-1 2 -1 0 0
0 -1 2 -1 0
0 0 -1 2 -1
0 0 0 -1 2

◆ getIndexer() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
__cuda_callable__ auto TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::getIndexer ( )

This method returns matrix elements indexer used by this matrix.

Returns
non-constant reference to the indexer.

◆ getIndexer() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
__cuda_callable__ auto TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::getIndexer ( ) const

This method returns matrix elements indexer used by this matrix.

Returns
constant reference to the indexer.

◆ getNonzeroElementsCount()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
Index TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::getNonzeroElementsCount ( ) const
overridevirtual

Returns number of non-zero matrix elements.

This method really counts the non-zero matrix elements and so it returns zero for matrix having all allocated elements set to zero.

Returns
number of non-zero matrix elements.

Reimplemented from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >.

◆ getRow() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
__cuda_callable__ auto TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::getRow ( IndexType rowIdx)

Non-constant getter of simple structure for accessing given matrix row.

Parameters
rowIdxis matrix row index.
Returns
RowView for accessing given matrix row.
Example
#include <iostream>
#include <TNL/Algorithms/parallelFor.h>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
getRowExample()
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 2 -1 . . . \ -> { 0, 0, 1 }
* | -1 2 -1 . . | -> { 0, 2, 1 }
* | . -1 2 -1. . | -> { 3, 2, 1 }
* | . . -1 2 -1 | -> { 3, 2, 1 }
* \ . . . -1 2 / -> { 3, 2, 1 }
*
*/
const int matrixSize( 5 );
MatrixType matrix( matrixSize, // number of matrix rows
matrixSize // number of matrix columns
);
auto view = matrix.getView();
auto f = [ = ] __cuda_callable__( int rowIdx ) mutable
{
auto row = view.getRow( rowIdx );
if( rowIdx > 0 )
row.setElement( 0, -1.0 ); // elements below the diagonal
row.setElement( 1, 2.0 ); // elements on the diagonal
if( rowIdx < matrixSize - 1 ) // elements above the diagonal
row.setElement( 2, -1.0 );
};
/***
* Set the matrix elements.
*/
TNL::Algorithms::parallelFor< Device >( 0, view.getRows(), f );
std::cout << std::endl << matrix << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Getting matrix rows on host: " << std::endl;
getRowExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Getting matrix rows on CUDA device: " << std::endl;
getRowExample< TNL::Devices::Cuda >();
#endif
}
Output
Getting matrix rows on host:
Row: 0 -> 0:2 1:-1
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 -> 3:-1 4:2
Getting matrix rows on CUDA device:
Row: 0 -> 0:2 1:-1
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 -> 3:-1 4:2

See TridiagonalMatrixRowView.

◆ getRow() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
__cuda_callable__ auto TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::getRow ( IndexType rowIdx) const

Constant getter of simple structure for accessing given matrix row.

Parameters
rowIdxis matrix row index.
Returns
RowView for accessing given matrix row.
Example
#include <iostream>
#include <functional>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
getRowExample()
{
const int matrixSize = 5;
MatrixType matrix( matrixSize, // number of matrix columns
{ { 0.0, 2.0, 1.0 }, // matrix elements
{ 0.0, 2.0, 1.0 },
{ 3.0, 2.0, 1.0 },
{ 3.0, 2.0, 1.0 },
{ 0.0, 2.0, 1.0 } } );
auto view = matrix.getView();
/***
* Fetch lambda function returns diagonal element in each row.
*/
auto fetch = [ = ] __cuda_callable__( int rowIdx ) mutable -> double
{
auto row = view.getRow( rowIdx );
return row.getValue( 2 ); // get value from subdiagonal with index 2, i.e. the main diagonal
};
/***
* Compute the matrix trace.
*/
int trace = TNL::Algorithms::reduce< Device >( 0, view.getRows(), fetch, std::plus<>{}, 0 );
std::cout << "Matrix reads as: " << std::endl << matrix << std::endl;
std::cout << "Matrix trace is: " << trace << "." << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Getting matrix rows on host: " << std::endl;
getRowExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Getting matrix rows on CUDA device: " << std::endl;
getRowExample< TNL::Devices::Cuda >();
#endif
}
Output
Getting matrix rows on host:
Matrix reads as:
Row: 0 -> 0:2 1:1
Row: 1 -> 1:2 2:1
Row: 2 -> 1:3 2:2 3:1
Row: 3 -> 2:3 3:2 4:1
Row: 4 -> 4:2
Matrix trace is: 5.
Getting matrix rows on CUDA device:
Matrix reads as:
Row: 0 -> 0:2 1:1
Row: 1 -> 1:2 2:1
Row: 2 -> 1:3 2:2 3:1
Row: 3 -> 2:3 3:2 4:1
Row: 4 -> 4:2
Matrix trace is: 5.

See TridiagonalMatrixRowView.

◆ getRowCapacities()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Vector >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::getRowCapacities ( Vector & rowCapacities) const

Compute capacities of all rows.

The row capacities are not stored explicitly and must be computed.

Parameters
rowCapacitiesis a vector where the row capacities will be stored.

◆ getSerializationType()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
std::string TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::getSerializationType ( )
static

Returns string with serialization type.

The string has a form Matrices::TridiagonalMatrix< RealType, [any_device], IndexType, Organization >.

See TridiagonalMatrix::getSerializationType.

Returns
String with the serialization type.

◆ operator!=()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_>
bool TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::operator!= ( const TridiagonalMatrixBase< Real_, Device_, Index_, Organization_ > & matrix) const

Comparison operator with another multidiagonal matrix.

Template Parameters
Real_is Real type of the source matrix.
Device_is Device type of the source matrix.
Index_is Index type of the source matrix.
Organization_is Organization of the source matrix.
Parameters
matrixis the source matrix view.
Returns
true if both matrices are NOT identical and false otherwise.

◆ operator=()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
TridiagonalMatrixBase & TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::operator= ( const TridiagonalMatrixBase< Real, Device, Index, Organization > & )
delete

Copy-assignment operator.

It is a deleted function, because matrix assignment in general requires reallocation.

◆ operator==()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization_>
bool TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::operator== ( const TridiagonalMatrixBase< Real_, Device_, Index_, Organization_ > & matrix) const

Comparison operator with another tridiagonal matrix.

Template Parameters
Real_is Real type of the source matrix.
Device_is Device type of the source matrix.
Index_is Index type of the source matrix.
Organization_is Organization of the source matrix.
Parameters
matrixis the source matrix view.
Returns
true if both matrices are identical and false otherwise.

◆ print()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::print ( std::ostream & str) const

Method for printing the matrix to output stream.

Parameters
stris the output stream.

◆ reduceAllRows()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::reduceAllRows ( Fetch & fetch,
Reduce & reduce,
Keep & keep,
const FetchReal & identity ) const

Method for performing general reduction on all matrix rows for constant instances.

Template Parameters
Fetchis a type of lambda function for data fetch declared as
auto fetch = [] __cuda_callable__ ( IndexType rowIdx, IndexType& columnIdx, RealType& elementValue ) -> FetchValue
{ ... };

The return type of this lambda can be any non void.

Template Parameters
Reduceis a type of lambda function for reduction declared as
auto reduce = [] __cuda_callable__ ( const FetchValue& v1, const FetchValue& v2 ) -> FetchValue { ... };
Template Parameters
Keepis a type of lambda function for storing results of reduction in each row. It is declared as
auto keep = [=] __cuda_callable__ ( IndexType rowIdx, const RealType& value ) { ... };
Template Parameters
FetchValueis type returned by the Fetch lambda function.
Parameters
fetchis an instance of lambda function for data fetch.
reduceis an instance of lambda function for reduction.
keepin an instance of lambda function for storing results.
identityis the identity element for the reduction operation, i.e. element which does not change the result of the reduction.
Example
#include <iostream>
#include <iomanip>
#include <functional>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void
reduceRows()
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 1 3 . . . \ -> { 0, 1, 3 }
* | 2 1 3 . . | -> { 2, 1, 3 }
* | . 2 1 3 . | -> { 2, 1, 3 }
* | . . 2 1 3 | -> { 2, 1, 3 }
* \ . . . 2 1 / 0 -> { 2, 1, 0 }
*
*/
TNL::Matrices::TridiagonalMatrix< double, Device > matrix( 5, // number of matrix columns
{ { 0, 1, 3 }, // matrix elements
{ 2, 1, 3 },
{ 2, 1, 3 },
{ 2, 1, 3 },
{ 2, 1, 3 } } );
auto view = matrix.getView();
/***
* Find largest element in each row.
*/
TNL::Containers::Vector< double, Device > rowMax( matrix.getRows() );
/***
* Prepare vector view for lambdas.
*/
auto rowMaxView = rowMax.getView();
/***
* Fetch lambda just returns absolute value of matrix elements.
*/
auto fetch = [] __cuda_callable__( int rowIdx, int columnIdx, const double& value ) -> double
{
return TNL::abs( value );
};
/***
* Reduce lambda return maximum of given values.
*/
auto reduce = [] __cuda_callable__( const double& a, const double& b ) -> double
{
return TNL::max( a, b );
};
/***
* Keep lambda store the largest value in each row to the vector rowMax.
*/
auto keep = [ = ] __cuda_callable__( int rowIdx, const double& value ) mutable
{
rowMaxView[ rowIdx ] = value;
};
/***
* Compute the largest values in each row.
*/
view.reduceAllRows( fetch, reduce, keep, std::numeric_limits< double >::lowest() ); // or matrix.reduceAllRows
std::cout << "The matrix reads as: " << std::endl << matrix << std::endl;
std::cout << "Max. elements in rows are: " << rowMax << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Rows reduction on host:" << std::endl;
reduceRows< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Rows reduction on CUDA device:" << std::endl;
reduceRows< TNL::Devices::Cuda >();
#endif
}
void reduceAllRows(Fetch &fetch, Reduce &reduce, Keep &keep, const FetchReal &identity) const
Method for performing general reduction on all matrix rows for constant instances.
Definition TridiagonalMatrixBase.hpp:229
Result reduce(Index begin, Index end, Fetch &&fetch, Reduction &&reduction, const Result &identity)
reduce implements (parallel) reduction for vectors and arrays.
Definition reduce.h:65
__cuda_callable__ T abs(const T &n)
This function returns absolute value of given number n.
Definition Math.h:74
constexpr ResultType max(const T1 &a, const T2 &b)
This function returns maximum of two numbers.
Definition Math.h:48
Output
Rows reduction on host:
The matrix reads as:
Row: 0 -> 0:1 1:3
Row: 1 -> 0:2 1:1 2:3
Row: 2 -> 1:2 2:1 3:3
Row: 3 -> 2:2 3:1 4:3
Row: 4 -> 3:2 4:1
Max. elements in rows are: [ 3, 3, 3, 3, 2 ]
Rows reduction on CUDA device:
The matrix reads as:
Row: 0 -> 0:1 1:3
Row: 1 -> 0:2 1:1 2:3
Row: 2 -> 1:2 2:1 3:3
Row: 3 -> 2:2 3:1 4:3
Row: 4 -> 3:2 4:1
Max. elements in rows are: [ 3, 3, 3, 3, 2 ]

◆ reduceRows()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::reduceRows ( IndexType begin,
IndexType end,
Fetch & fetch,
Reduce & reduce,
Keep & keep,
const FetchReal & identity ) const

Method for performing general reduction on matrix rows for constant instances.

Template Parameters
Fetchis a type of lambda function for data fetch declared as
auto fetch = [] __cuda_callable__ ( IndexType rowIdx, IndexType& columnIdx, RealType& elementValue ) -> FetchValue
{ ... };

The return type of this lambda can be any non void.

Template Parameters
Reduceis a type of lambda function for reduction declared as
auto reduce = [] __cuda_callable__ ( const FetchValue& v1, const FetchValue& v2 ) -> FetchValue { ... };
Template Parameters
Keepis a type of lambda function for storing results of reduction in each row. It is declared as
auto keep = [=] __cuda_callable__ ( IndexType rowIdx, const RealType& value ) { ... };
Template Parameters
FetchValueis type returned by the Fetch lambda function.
Parameters
begindefines beginning of the range [begin, end) of rows to be processed.
enddefines ending of the range [begin, end) of rows to be processed.
fetchis an instance of lambda function for data fetch.
reduceis an instance of lambda function for reduction.
keepin an instance of lambda function for storing results.
identityis the identity element for the reduction operation, i.e. element which does not change the result of the reduction.
Example
#include <iostream>
#include <iomanip>
#include <functional>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void
reduceRows()
{
/***
* Set the following matrix (dots represent zero matrix elements and zeros are
* padding zeros for memory alignment):
*
* 0 / 1 3 . . . \ -> { 0, 1, 3 }
* | 2 1 3 . . | -> { 2, 1, 3 }
* | . 2 1 3 . | -> { 2, 1, 3 }
* | . . 2 1 3 | -> { 2, 1, 3 }
* \ . . . 2 1 / 0 -> { 2, 1, 0 }
*
*/
TNL::Matrices::TridiagonalMatrix< double, Device > matrix( 5, // number of matrix columns
{ { 0, 1, 3 }, // matrix elements
{ 2, 1, 3 },
{ 2, 1, 3 },
{ 2, 1, 3 },
{ 2, 1, 3 } } );
auto view = matrix.getView();
/***
* Find largest element in each row.
*/
TNL::Containers::Vector< double, Device > rowMax( matrix.getRows() );
/***
* Prepare vector view for lambdas.
*/
auto rowMaxView = rowMax.getView();
/***
* Fetch lambda just returns absolute value of matrix elements.
*/
auto fetch = [] __cuda_callable__( int rowIdx, int columnIdx, const double& value ) -> double
{
return TNL::abs( value );
};
/***
* Reduce lambda return maximum of given values.
*/
auto reduce = [] __cuda_callable__( const double& a, const double& b ) -> double
{
return TNL::max( a, b );
};
/***
* Keep lambda store the largest value in each row to the vector rowMax.
*/
auto keep = [ = ] __cuda_callable__( int rowIdx, const double& value ) mutable
{
rowMaxView[ rowIdx ] = value;
};
/***
* Compute the largest values in each row.
*/
view.reduceRows( 0, view.getRows(), fetch, reduce, keep, std::numeric_limits< double >::lowest() ); // or matrix.reduceRows
std::cout << "The matrix reads as: " << std::endl << matrix << std::endl;
std::cout << "Max. elements in rows are: " << rowMax << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Rows reduction on host:" << std::endl;
reduceRows< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Rows reduction on CUDA device:" << std::endl;
reduceRows< TNL::Devices::Cuda >();
#endif
}
__cuda_callable__ IndexType getRows() const
Returns number of matrix rows.
Definition MatrixBase.hpp:52
void reduceRows(IndexType begin, IndexType end, Fetch &fetch, Reduce &reduce, Keep &keep, const FetchReal &identity) const
Method for performing general reduction on matrix rows for constant instances.
Definition TridiagonalMatrixBase.hpp:188
Output
Rows reduction on host:
The matrix reads as:
Row: 0 -> 0:1 1:3
Row: 1 -> 0:2 1:1 2:3
Row: 2 -> 1:2 2:1 3:3
Row: 3 -> 2:2 3:1 4:3
Row: 4 -> 3:2 4:1
Max. elements in rows are: [ 3, 3, 3, 3, 2 ]
Rows reduction on CUDA device:
The matrix reads as:
Row: 0 -> 0:1 1:3
Row: 1 -> 0:2 1:1 2:3
Row: 2 -> 1:2 2:1 3:3
Row: 3 -> 2:2 3:1 4:3
Row: 4 -> 3:2 4:1
Max. elements in rows are: [ 3, 3, 3, 3, 2 ]

◆ sequentialForAllRows() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Function >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::sequentialForAllRows ( Function & function)

This method calls sequentialForRows for all matrix rows.

See TridiagonalMatrixBase::sequentialForAllRows.

Template Parameters
Functionis a type of lambda function that will operate on matrix elements.
Parameters
functionis an instance of the lambda function to be called in each row.

◆ sequentialForAllRows() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Function >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::sequentialForAllRows ( Function & function) const

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

See TridiagonalMatrixBase::sequentialForRows.

Template Parameters
Functionis a type of lambda function that will operate on matrix elements.
Parameters
functionis an instance of the lambda function to be called in each row.

◆ sequentialForRows() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Function >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::sequentialForRows ( IndexType begin,
IndexType end,
Function & function )

Method for sequential iteration over all matrix rows for non-constant instances.

Template Parameters
Functionis type of lambda function that will operate on matrix elements. It should have form like
auto function = [] __cuda_callable__ ( RowView& row ) { ... };

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

Parameters
begindefines beginning of the range [begin, end) of rows to be processed.
enddefines ending of the range [begin, end) of rows to be processed.
functionis an instance of the lambda function to be called in each row.

◆ sequentialForRows() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename Function >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::sequentialForRows ( IndexType begin,
IndexType end,
Function & function ) const

Method for sequential iteration over all matrix rows for constant instances.

Template Parameters
Functionis type of lambda function that will operate on matrix elements. It should have form like
auto function = [] __cuda_callable__ ( const ConstRowView& row ) { ... };

ConstRowView represents matrix row - see TNL::Matrices::TridiagonalMatrixBase::ConstRowView.

Parameters
begindefines beginning of the range [begin, end) of rows to be processed.
enddefines ending of the range [begin, end) of rows to be processed.
functionis an instance of the lambda function to be called in each row.

◆ setElement()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
__cuda_callable__ void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::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. TridiagonalMatrix::getRow or TridiagonalMatrix::forElements and TridiagonalMatrix::forAllElements. The call may fail if the matrix row capacity is exhausted.

Parameters
rowis row index of the element.
columnis columns index of the element.
valueis the value the element will be set to.
Example
#include <iostream>
#include <TNL/Algorithms/parallelFor.h>
#include <TNL/Matrices/TridiagonalMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void
setElements()
{
const int matrixSize( 5 );
Matrix matrix( matrixSize, matrixSize );
auto view = matrix.getView();
for( int i = 0; i < 5; i++ )
view.setElement( i, i, i ); // or matrix.setElement
std::cout << "Matrix set from the host:" << std::endl;
std::cout << matrix << std::endl;
auto f = [ = ] __cuda_callable__( int i ) mutable
{
if( i > 0 )
view.setElement( i, i - 1, 1.0 );
view.setElement( i, i, -i );
if( i < matrixSize - 1 )
view.setElement( i, i + 1, 1.0 );
};
TNL::Algorithms::parallelFor< Device >( 0, matrixSize, f );
std::cout << "Matrix set from its native device:" << std::endl;
std::cout << matrix << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Set elements on host:" << std::endl;
setElements< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Set elements on CUDA device:" << std::endl;
setElements< TNL::Devices::Cuda >();
#endif
}
Output
Set elements on host:
Matrix set from the host:
Row: 0 ->
Row: 1 -> 1:1
Row: 2 -> 2:2
Row: 3 -> 3:3
Row: 4 -> 4:4
Matrix set from its native device:
Row: 0 -> 1:1
Row: 1 -> 0:1 1:-1 2:1
Row: 2 -> 1:1 2:-2 3:1
Row: 3 -> 2:1 3:-3 4:1
Row: 4 -> 3:1 4:-4
Set elements on CUDA device:
Matrix set from the host:
Row: 0 ->
Row: 1 -> 1:1
Row: 2 -> 2:2
Row: 3 -> 3:3
Row: 4 -> 4:4
Matrix set from its native device:
Row: 0 -> 1:1
Row: 1 -> 0:1 1:-1 2:1
Row: 2 -> 1:1 2:-2 3:1
Row: 3 -> 2:1 3:-3 4:1
Row: 4 -> 3:1 4:-4

◆ setValue()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::setValue ( const RealType & value)

Set all matrix elements to given value.

Parameters
valueis the new value of all matrix elements.

◆ vectorProduct()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization>
template<typename InVector , typename OutVector >
void TNL::Matrices::TridiagonalMatrixBase< Real, Device, Index, Organization >::vectorProduct ( const InVector & inVector,
OutVector & outVector,
RealType matrixMultiplicator = 1.0,
RealType outVectorMultiplicator = 0.0,
IndexType begin = 0,
IndexType end = 0 ) const

Computes product of matrix and vector.

More precisely, it computes:

outVector = matrixMultiplicator * ( * this ) * inVector + outVectorMultiplicator * outVector
Template Parameters
InVectoris type of input vector. It can be TNL::Containers::Vector, TNL::Containers::VectorView, TNL::Containers::Array, TNL::Containers::ArrayView, or similar container.
OutVectoris type of output vector. It can be TNL::Containers::Vector, TNL::Containers::VectorView, TNL::Containers::Array, TNL::Containers::ArrayView, or similar container.
Parameters
inVectoris input vector.
outVectoris output vector.
matrixMultiplicatoris a factor by which the matrix is multiplied. It is one by default.
outVectorMultiplicatoris a factor by which the outVector is multiplied before added to the result of matrix-vector product. It is zero by default.
beginis the beginning of the rows range for which the vector product is computed. It is zero by default.
endis the end of the rows range for which the vector product is computed. It is number if the matrix rows by default.

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