Template Numerical Library version\ main:9e7b0f4
Loading...
Searching...
No Matches
TNL::Solvers::Linear::UmfpackWrapper< Matrix, SolverMonitor > Class Template Reference
Inheritance diagram for TNL::Solvers::Linear::UmfpackWrapper< Matrix, SolverMonitor >:
Collaboration diagram for TNL::Solvers::Linear::UmfpackWrapper< Matrix, SolverMonitor >:

Public Types

using ConstVectorViewType = typename Base::ConstVectorViewType
using DeviceType = typename Base::DeviceType
using IndexType = typename Base::IndexType
using MatrixPointer = typename Base::MatrixPointer
using MatrixType = typename Base::MatrixType
using RealType = typename Base::RealType
using VectorViewType = typename Base::VectorViewType
Public Types inherited from TNL::Solvers::Linear::LinearSolver< Matrix >
using ConstVectorViewType = typename Traits< Matrix >::ConstVectorViewType
 Type for constant vector view.
using DeviceType = typename Matrix::DeviceType
 Device where the solver will run on and auxiliary 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

void setMatrix (const MatrixPointer &matrix) override
 Set the matrix of the linear system.
bool solve (ConstVectorViewType b, VectorViewType x) override
 Method for solving a linear system.
Public Member Functions inherited from TNL::Solvers::Linear::LinearSolver< Matrix >
virtual ~LinearSolver ()=default
 Default destructor.
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.
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.

Protected Attributes

bool factorized = false
MatrixPointer matrix
Protected Attributes inherited from TNL::Solvers::Linear::LinearSolver< Matrix >
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< double > * solverMonitor

Additional Inherited Members

Static Public Member Functions inherited from TNL::Solvers::Linear::LinearSolver< Matrix >
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.

Member Function Documentation

◆ setMatrix()

template<typename Matrix, typename SolverMonitor>
void TNL::Solvers::Linear::UmfpackWrapper< Matrix, SolverMonitor >::setMatrix ( const MatrixPointer & matrix)
overridevirtual

Set the matrix of the linear system.

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

Reimplemented from TNL::Solvers::Linear::LinearSolver< Matrix >.

◆ solve()

template<typename Matrix, typename SolverMonitor>
bool TNL::Solvers::Linear::UmfpackWrapper< Matrix, SolverMonitor >::solve ( ConstVectorViewType b,
VectorViewType x )
overridevirtual

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
}
#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:55
virtual void setMatrix(const MatrixPointer &matrix)
Set the matrix of the linear system.
Definition LinearSolver.h:120
Iterative solver of linear systems based on the Transpose-free quasi-minimal residual (TFQMR) method.
Definition TFQMR.h:21
T endl(T... args)
T make_shared(T... args)
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 ]

Implements TNL::Solvers::Linear::LinearSolver< Matrix >.


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