Template Numerical Library version\ main:9e7b0f4
Loading...
Searching...
No Matches
TNL::Solvers::Linear::GEM< Matrix, Real, SolverMonitor > Struct Template Reference

Gaussian Elimination Method (GEM) direct solver for dense linear systems. More...

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

Inheritance diagram for TNL::Solvers::Linear::GEM< Matrix, Real, SolverMonitor >:
Collaboration diagram for TNL::Solvers::Linear::GEM< Matrix, Real, 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

 GEM ()=default
 Default constructor.
bool getPivoting () const
 Checks if pivoting is enabled.
void setPivoting (bool pivoting)
 Enables or disables pivoting in the GEM algorithm.
bool solve (ConstVectorViewType b, VectorViewType x) override
 Solves the linear system Ax = b.
bool solve (MatrixType &A, ConstVectorViewType b, VectorViewType x)
 Solves the linear system Ax = b.
Public Member Functions inherited from TNL::Solvers::Linear::LinearSolver< Matrix >
virtual ~LinearSolver ()=default
 Default destructor.
virtual 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.
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 pivoting = true
 Indicates whether pivoting is enabled.
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.

Detailed Description

template<typename Matrix, typename Real = typename Matrix::RealType, typename SolverMonitor = IterativeSolverMonitor< double >>
struct TNL::Solvers::Linear::GEM< Matrix, Real, SolverMonitor >

Gaussian Elimination Method (GEM) direct solver for dense linear systems.

This class implements the GEM algorithm for solving linear systems with dense matrices. It supports optional pivoting.

Template Parameters
MatrixType of the matrix representing the linear system. Must be a dense matrix.
RealFloating point type used for computations.
SolverMonitorType of the solver monitor.

Member Function Documentation

◆ getPivoting()

template<typename Matrix, typename Real, typename SolverMonitor>
bool TNL::Solvers::Linear::GEM< Matrix, Real, SolverMonitor >::getPivoting ( ) const
nodiscard

Checks if pivoting is enabled.

Returns
True if pivoting is enabled, false otherwise.

◆ setPivoting()

template<typename Matrix, typename Real, typename SolverMonitor>
void TNL::Solvers::Linear::GEM< Matrix, Real, SolverMonitor >::setPivoting ( bool pivoting)

Enables or disables pivoting in the GEM algorithm.

Parameters
pivotingIf true, pivoting is enabled; otherwise, it is disabled.

◆ solve() [1/2]

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

Solves the linear system Ax = b.

This method makes a copy of the matrix A. In case that it is consuming too much memory use method solve taking matrix as input parameter.

Parameters
bRight-hand side vector.
xSolution vector (output).
Returns
True if the solver succeeded, false otherwise.

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

◆ solve() [2/2]

template<typename Matrix, typename Real, typename SolverMonitor>
bool TNL::Solvers::Linear::GEM< Matrix, Real, SolverMonitor >::solve ( MatrixType & A,
ConstVectorViewType b,
VectorViewType x )

Solves the linear system Ax = b.

Parameters
AMatrix representing the linear system. This matrix is modified during the execution of this method.
bRight-hand side vector.
xSolution vector (output).
Returns
True if the solver succeeded, false otherwise.

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