Template Numerical Library version\ main:1437bf49
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Protected Attributes | List of all members
TNL::Solvers::Linear::LinearSolver< Matrix > Class Template Referenceabstract

Base class for iterative solvers of systems of linear equations. More...

#include <TNL/Solvers/Linear/LinearSolver.h>

Inheritance diagram for TNL::Solvers::Linear::LinearSolver< Matrix >:
Inheritance graph
[legend]
Collaboration diagram for TNL::Solvers::Linear::LinearSolver< Matrix >:
Collaboration graph
[legend]

Public Types

using ConstVectorViewType = typename Traits< Matrix >::ConstVectorViewType
 Type for constant vector view.
 
using DeviceType = typename Matrix::DeviceType
 Device where the solver will run on and auxillary data will alloacted on.
 
using IndexType = typename Matrix::IndexType
 Type for indexing.
 
using MatrixPointer = std::shared_ptr< std::add_const_t< MatrixType > >
 Type of shared pointer to the matrix.
 
using MatrixType = Matrix
 Type of the matrix representing the linear system.
 
using PreconditionerPointer = std::shared_ptr< std::add_const_t< PreconditionerType > >
 Type of shared pointer to the preconditioner.
 
using PreconditionerType = Preconditioners::Preconditioner< MatrixType >
 Type of preconditioner.
 
using RealType = typename Matrix::RealType
 Floating point type used for computations.
 
using VectorViewType = typename Traits< Matrix >::VectorViewType
 Type for vector view.
 
- Public Types inherited from TNL::Solvers::IterativeSolver< Matrix::RealType, Matrix::IndexType >
using SolverMonitorType
 Type of an object used for monitoring of the convergence.
 

Public Member Functions

virtual ~LinearSolver ()
 Default destructor.
 
void setMatrix (const MatrixPointer &matrix)
 Set the matrix of the linear system.
 
void setPreconditioner (const PreconditionerPointer &preconditioner)
 Set the preconditioner.
 
virtual bool setup (const Config::ParameterContainer &parameters, const String &prefix="")
 Method for setup of the linear iterative solver based on configuration parameters.
 
virtual bool solve (ConstVectorViewType b, VectorViewType x)=0
 Method for solving a linear system.
 
- Public Member Functions inherited from TNL::Solvers::IterativeSolver< Matrix::RealType, Matrix::IndexType >
 IterativeSolver ()=default
 Default constructor.
 
bool checkConvergence ()
 Checks whether the convergence occurred already.
 
bool checkNextIteration ()
 Checks if the solver is allowed to do the next iteration.
 
const Matrix::RealType & getConvergenceResidue () const
 Gets the the convergence threshold.
 
const Matrix::RealType & getDivergenceResidue () const
 Gets the limit for the divergence criterion.
 
const Matrix::IndexType & getIterations () const
 Gets the number of iterations performed by the solver so far.
 
const Matrix::IndexType & getMaxIterations () const
 Gets the maximal number of iterations the solver is allowed to perform.
 
const Matrix::IndexType & getMinIterations () const
 Gets the minimal number of iterations the solver is supposed to do.
 
const Matrix::RealType & getResidue () const
 Gets the residue reached at the current iteration.
 
bool nextIteration ()
 Proceeds to the next iteration.
 
void resetIterations ()
 Sets the the number of the current iterations to zero.
 
void setConvergenceResidue (const Matrix::RealType &convergenceResidue)
 Sets the threshold for the convergence.
 
void setDivergenceResidue (const Matrix::RealType &divergenceResidue)
 Sets the residue limit for the divergence criterion.
 
void setMaxIterations (const Matrix::IndexType &maxIterations)
 Sets the maximal number of iterations the solver is allowed to perform.
 
void setMinIterations (const Matrix::IndexType &minIterations)
 Sets the minimal number of iterations the solver is supposed to do.
 
void setRefreshRate (const Matrix::IndexType &refreshRate)
 Sets the refresh rate (in milliseconds) for the solver monitor.
 
void setResidue (const Matrix::RealType &residue)
 Sets the residue reached at the current iteration.
 
void setSolverMonitor (SolverMonitorType &solverMonitor)
 Sets the solver monitor object.
 
bool setup (const Config::ParameterContainer &parameters, const std::string &prefix="")
 Method for setup of the iterative solver based on configuration parameters.
 

Static Public Member Functions

static void configSetup (Config::ConfigDescription &config, const String &prefix="")
 This method defines configuration entries for setup of the linear iterative solver.
 
- Static Public Member Functions inherited from TNL::Solvers::IterativeSolver< Matrix::RealType, Matrix::IndexType >
static void configSetup (Config::ConfigDescription &config, const std::string &prefix="")
 This method defines configuration entries for setup of the iterative solver.
 

Protected Attributes

MatrixPointer matrix = nullptr
 
PreconditionerPointer preconditioner = nullptr
 
- Protected Attributes inherited from TNL::Solvers::IterativeSolver< Matrix::RealType, Matrix::IndexType >
Matrix::RealType convergenceResidue
 
Matrix::IndexType currentIteration
 
Matrix::RealType currentResidue
 
Matrix::RealType divergenceResidue
 
Matrix::IndexType maxIterations
 
Matrix::IndexType minIterations
 
Matrix::IndexType refreshRate
 
std::ofstream residualHistoryFile
 
std::string residualHistoryFileName
 
IterativeSolverMonitor< Matrix::RealType, Matrix::IndexType > * solverMonitor
 

Detailed Description

template<typename Matrix>
class TNL::Solvers::Linear::LinearSolver< Matrix >

Base class for iterative solvers of systems of linear equations.

To use the linear solver, one needs to first set the matrix of the linear system by means of the method LinearSolver::setMatrix. Afterward, one may call the method LinearSolver::solve which accepts the right-hand side vector b and a vector x to which the solution will be stored. One may also use appropriate preconditioner to speed-up the convergence - see the method LinearSolver::setPreconditioner.

Template Parameters
Matrixis type of matrix representing the linear system.

The following example demonstrates the use iterative linear solvers:

1#include <iostream>
2#include <memory>
3#include <TNL/Matrices/SparseMatrix.h>
4#include <TNL/Devices/Host.h>
5#include <TNL/Devices/Cuda.h>
6#include <TNL/Solvers/Linear/TFQMR.h>
7
8template< typename Device >
9void
10iterativeLinearSolverExample()
11{
12 /***
13 * Set the following matrix (dots represent zero matrix elements):
14 *
15 * / 2.5 -1 . . . \
16 * | -1 2.5 -1 . . |
17 * | . -1 2.5 -1. . |
18 * | . . -1 2.5 -1 |
19 * \ . . . -1 2.5 /
20 */
23 const int size( 5 );
24 auto matrix_ptr = std::make_shared< MatrixType >();
25 matrix_ptr->setDimensions( size, size );
26 matrix_ptr->setRowCapacities( Vector( { 2, 3, 3, 3, 2 } ) );
27
28 auto f = [ = ] __cuda_callable__( typename MatrixType::RowView & row ) mutable
29 {
30 const int rowIdx = row.getRowIndex();
31 if( rowIdx == 0 ) {
32 row.setElement( 0, rowIdx, 2.5 ); // diagonal element
33 row.setElement( 1, rowIdx + 1, -1 ); // element above the diagonal
34 }
35 else if( rowIdx == size - 1 ) {
36 row.setElement( 0, rowIdx - 1, -1.0 ); // element below the diagonal
37 row.setElement( 1, rowIdx, 2.5 ); // diagonal element
38 }
39 else {
40 row.setElement( 0, rowIdx - 1, -1.0 ); // element below the diagonal
41 row.setElement( 1, rowIdx, 2.5 ); // diagonal element
42 row.setElement( 2, rowIdx + 1, -1.0 ); // element above the diagonal
43 }
44 };
45
46 /***
47 * Set the matrix elements.
48 */
49 matrix_ptr->forAllRows( f );
50 std::cout << *matrix_ptr << std::endl;
51
52 /***
53 * Set the right-hand side vector.
54 */
55 Vector x( size, 1.0 );
56 Vector b( size );
57 matrix_ptr->vectorProduct( x, b );
58 x = 0.0;
59 std::cout << "Vector b = " << b << std::endl;
60
61 /***
62 * Solve the linear system.
63 */
65 LinearSolver solver;
66 solver.setMatrix( matrix_ptr );
67 solver.setConvergenceResidue( 1.0e-6 );
68 solver.solve( b, x );
69 std::cout << "Vector x = " << x << std::endl;
70}
71
72int
73main( int argc, char* argv[] )
74{
75 std::cout << "Solving linear system on host: " << std::endl;
76 iterativeLinearSolverExample< TNL::Devices::Sequential >();
77
78#ifdef __CUDACC__
79 std::cout << "Solving linear system on CUDA device: " << std::endl;
80 iterativeLinearSolverExample< TNL::Devices::Cuda >();
81#endif
82}
#define __cuda_callable__
Definition Macros.h:49
Vector extends Array with algebraic operations.
Definition Vector.h:36
Implementation of sparse matrix, i.e. matrix storing only non-zero elements.
Definition SparseMatrix.h:57
void setMatrix(const MatrixPointer &matrix)
Set the matrix of the linear system.
Definition LinearSolver.h:120
Matrix MatrixType
Type of the matrix representing the linear system.
Definition LinearSolver.h:71
Iterative solver of linear systems based on the Transpose-free quasi-minimal residual (TFQMR) method.
Definition TFQMR.h:21
T endl(T... args)

The result looks as follows:

Solving linear system on host:
Row: 0 -> 0:2.5 1:-1
Row: 1 -> 0:-1 1:2.5 2:-1
Row: 2 -> 1:-1 2:2.5 3:-1
Row: 3 -> 2:-1 3:2.5 4:-1
Row: 4 -> 3:-1 4:2.5
Vector b = [ 1.5, 0.5, 0.5, 0.5, 1.5 ]
Vector x = [ 1, 1, 1, 1, 1 ]
Solving linear system on CUDA device:
Row: 0 -> 0:2.5 1:-1
Row: 1 -> 0:-1 1:2.5 2:-1
Row: 2 -> 1:-1 2:2.5 3:-1
Row: 3 -> 2:-1 3:2.5 4:-1
Row: 4 -> 3:-1 4:2.5
Vector b = [ 1.5, 0.5, 0.5, 0.5, 1.5 ]
Vector x = [ 1, 1, 1, 1, 1 ]

See also TNL::Solvers::IterativeSolverMonitor for monitoring of iterative solvers.

Member Typedef Documentation

◆ DeviceType

template<typename Matrix >
using TNL::Solvers::Linear::LinearSolver< Matrix >::DeviceType = typename Matrix::DeviceType

Device where the solver will run on and auxillary data will alloacted on.

See Devices::Host or Devices::Cuda.

Member Function Documentation

◆ configSetup()

template<typename Matrix >
static void TNL::Solvers::Linear::LinearSolver< Matrix >::configSetup ( Config::ConfigDescription & config,
const String & prefix = "" )
inlinestatic

This method defines configuration entries for setup of the linear iterative solver.

See IterativeSolver::configSetup.

Parameters
configcontains description of configuration parameters.
prefixis a prefix of particular configuration entries.

◆ setMatrix()

template<typename Matrix >
void TNL::Solvers::Linear::LinearSolver< Matrix >::setMatrix ( const MatrixPointer & matrix)
inline

Set the matrix of the linear system.

Parameters
matrixis a shared pointer to the matrix of the linear system

◆ setPreconditioner()

template<typename Matrix >
void TNL::Solvers::Linear::LinearSolver< Matrix >::setPreconditioner ( const PreconditionerPointer & preconditioner)
inline

Set the preconditioner.

Parameters
preconditioneris a shared pointer to preconditioner.

◆ setup()

template<typename Matrix >
virtual bool TNL::Solvers::Linear::LinearSolver< Matrix >::setup ( const Config::ParameterContainer & parameters,
const String & prefix = "" )
inlinevirtual

Method for setup of the linear iterative solver based on configuration parameters.

Parameters
parameterscontains values of the define configuration entries.
prefixis a prefix of particular configuration entries.

Reimplemented in TNL::Solvers::Linear::BICGStab< Matrix >, TNL::Solvers::Linear::BICGStabL< Matrix >, TNL::Solvers::Linear::GMRES< Matrix >, TNL::Solvers::Linear::IDRs< Matrix >, TNL::Solvers::Linear::Jacobi< Matrix >, and TNL::Solvers::Linear::SOR< Matrix >.

◆ solve()

template<typename Matrix >
virtual bool TNL::Solvers::Linear::LinearSolver< Matrix >::solve ( ConstVectorViewType b,
VectorViewType x )
pure virtual

Method for solving a linear system.

The linear system is defined by the matrix given by the method LinearSolver::setMatrix and by the right-hand side vector represented by the vector b. The result is stored in the vector b. The solver can be accelerated with appropriate preconditioner set by the methods LinearSolver::setPreconditioner.

Parameters
bvector with the right-hand side of the linear system.
xvector for the solution of the linear system.
Returns
true if the solver converged.
false if the solver did not converge.
Example
#include <iostream>
#include <memory>
#include <TNL/Matrices/SparseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
#include <TNL/Solvers/Linear/TFQMR.h>
template< typename Device >
void
iterativeLinearSolverExample()
{
/***
* Set the following matrix (dots represent zero matrix elements):
*
* / 2.5 -1 . . . \
* | -1 2.5 -1 . . |
* | . -1 2.5 -1. . |
* | . . -1 2.5 -1 |
* \ . . . -1 2.5 /
*/
const int size( 5 );
auto matrix_ptr = std::make_shared< MatrixType >();
matrix_ptr->setDimensions( size, size );
matrix_ptr->setRowCapacities( Vector( { 2, 3, 3, 3, 2 } ) );
auto f = [ = ] __cuda_callable__( typename MatrixType::RowView & row ) mutable
{
const int rowIdx = row.getRowIndex();
if( rowIdx == 0 ) {
row.setElement( 0, rowIdx, 2.5 ); // diagonal element
row.setElement( 1, rowIdx + 1, -1 ); // element above the diagonal
}
else if( rowIdx == size - 1 ) {
row.setElement( 0, rowIdx - 1, -1.0 ); // element below the diagonal
row.setElement( 1, rowIdx, 2.5 ); // diagonal element
}
else {
row.setElement( 0, rowIdx - 1, -1.0 ); // element below the diagonal
row.setElement( 1, rowIdx, 2.5 ); // diagonal element
row.setElement( 2, rowIdx + 1, -1.0 ); // element above the diagonal
}
};
/***
* Set the matrix elements.
*/
matrix_ptr->forAllRows( f );
std::cout << *matrix_ptr << std::endl;
/***
* Set the right-hand side vector.
*/
Vector x( size, 1.0 );
Vector b( size );
matrix_ptr->vectorProduct( x, b );
x = 0.0;
std::cout << "Vector b = " << b << std::endl;
/***
* Solve the linear system.
*/
LinearSolver solver;
solver.setMatrix( matrix_ptr );
solver.setConvergenceResidue( 1.0e-6 );
solver.solve( b, x );
std::cout << "Vector x = " << x << std::endl;
}
int
main( int argc, char* argv[] )
{
std::cout << "Solving linear system on host: " << std::endl;
iterativeLinearSolverExample< TNL::Devices::Sequential >();
#ifdef __CUDACC__
std::cout << "Solving linear system on CUDA device: " << std::endl;
iterativeLinearSolverExample< TNL::Devices::Cuda >();
#endif
}
Output
Solving linear system on host:
Row: 0 -> 0:2.5 1:-1
Row: 1 -> 0:-1 1:2.5 2:-1
Row: 2 -> 1:-1 2:2.5 3:-1
Row: 3 -> 2:-1 3:2.5 4:-1
Row: 4 -> 3:-1 4:2.5
Vector b = [ 1.5, 0.5, 0.5, 0.5, 1.5 ]
Vector x = [ 1, 1, 1, 1, 1 ]
Solving linear system on CUDA device:
Row: 0 -> 0:2.5 1:-1
Row: 1 -> 0:-1 1:2.5 2:-1
Row: 2 -> 1:-1 2:2.5 3:-1
Row: 3 -> 2:-1 3:2.5 4:-1
Row: 4 -> 3:-1 4:2.5
Vector b = [ 1.5, 0.5, 0.5, 0.5, 1.5 ]
Vector x = [ 1, 1, 1, 1, 1 ]

Implemented in TNL::Solvers::Linear::BICGStab< Matrix >, TNL::Solvers::Linear::BICGStabL< Matrix >, TNL::Solvers::Linear::CG< Matrix >, TNL::Solvers::Linear::GMRES< Matrix >, TNL::Solvers::Linear::IDRs< Matrix >, TNL::Solvers::Linear::Jacobi< Matrix >, TNL::Solvers::Linear::SOR< Matrix >, and TNL::Solvers::Linear::TFQMR< Matrix >.


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