NeoPZ
Classes | Macros | Functions
The Solver classes

This module contains all classes that represent a matrix inversion procedure.
Representing a inversion procedure as an object gives the user great flexibility to combine solution procedures. More...

Classes

class  TPZSequenceSolver< TVar >
 Defines sequence solvers. Solver. More...
 
class  TPZSolver< TVar >
 Defines a abstract class of solvers which will be used by matrix classes. Solver. More...
 
class  TPZMatrixSolver< TVar >
 Defines a class of matrix solvers. Solver. More...
 
class  TPZStepSolver< TVar >
 Defines step solvers class. Solver. More...
 
class  TPZCopySolve< TVar >
 To solve clones of the given matrix. Solver. More...
 
class  TPZMGSolver< TVar >
 Represents a solution process in three steps: transfer of the residual, execute a solver on the coarse mesh, extend the solution. Solver. More...
 

Macros

#define TPZSQUENCESOLVER_ID
 Id for sequence solver. More...
 
#define TPZMATRIXSOLVER_ID
 

Functions

template<class Matrix , class Vector , class Preconditioner , class Real >
int CG (Matrix &A, Vector &x, const Vector &b, Preconditioner &M, Vector *residual, int64_t &max_iter, Real &tol, const int FromCurrent)
 CG solves the symmetric positive definite linear system $ A x = b $ using the Conjugate Gradient method. More...
 
template<class Matrix , class Vector , class Preconditioner , class Real >
int CG (const Matrix &A, Vector &x, const Vector &b, const Preconditioner &M, int64_t &max_iter, Real &tol)
 CG solves the symmetric positive definite linear system $ Ax=b $ using the Conjugate Gradient method. More...
 
template<class Matrix , class Vector , class Preconditioner , class Real >
int CGS (const Matrix &A, Vector &x, const Vector &b, const Preconditioner &M, int64_t &max_iter, Real &tol)
 CGS solves the unsymmetric linear system $ Ax = b $ using the Conjugate Gradient Squared method. More...
 
template<class Matrix , class Vector , class Preconditioner , class Real , class Type >
int CHEBY (const Matrix &A, Vector &x, const Vector &b, const Preconditioner &M, int64_t &max_iter, Real &tol, Type eigmin, Type eigmax)
 CHEBY solves the symmetric positive definite linear system $ Ax = b $ using the Preconditioned Chebyshev Method. More...
 
template<class Matrix , class Vector >
void Update (Vector &x, int64_t k, Matrix &h, Vector &s, Vector v[])
 Computes back solve. Updates on v vector. More...
 
template<class Real >
Real abs (Real x)
 Returns absolute value to real. More...
 
template<class Operator , class Vector , class Preconditioner , class Matrix , class Real >
int GMRES (Operator &A, Vector &x, const Vector &b, Preconditioner &M, Matrix &H, int &m, int64_t &max_iter, Real &tol, Vector *residual, const int FromCurrent)
 GMRES solves the unsymmetric linear system Ax = b using the Generalized Minimum Residual method. More...
 
template<class Matrix , class Vector , class Preconditioner1 , class Preconditioner2 , class Real >
int QMR (const Matrix &A, Vector &x, const Vector &b, const Preconditioner1 &M1, const Preconditioner2 &M2, int64_t &max_iter, Real &tol)
 QMR solves the unsymmetric linear system $ Ax = b $ using the Quasi-Minimal Residual method. More...
 

Detailed Description

This module contains all classes that represent a matrix inversion procedure.
Representing a inversion procedure as an object gives the user great flexibility to combine solution procedures.

Macro Definition Documentation

◆ TPZMATRIXSOLVER_ID

#define TPZMATRIXSOLVER_ID

Definition at line 66 of file pzsolve.h.

◆ TPZSQUENCESOLVER_ID

#define TPZSQUENCESOLVER_ID

Id for sequence solver.

Definition at line 16 of file pzseqsolver.h.

Function Documentation

◆ abs()

template<class Real >
Real abs ( Real  x)

Returns absolute value to real.

Definition at line 43 of file gmres.h.

Referenced by GeneratePlaneRotation().

◆ CG() [1/2]

template<class Matrix , class Vector , class Preconditioner , class Real >
int CG ( const Matrix &  A,
Vector x,
const Vector b,
const Preconditioner &  M,
int64_t &  max_iter,
Real &  tol 
)

CG solves the symmetric positive definite linear system $ Ax=b $ using the Conjugate Gradient method.

Returns
The return value indicates convergence within max_iter (input) iterations (0), or no convergence within max_iter iterations (1).
Upon successful return, output arguments have the following values:
Parameters
A– matrix of the system
b– vector of the system
M– preconditioner matrix
x– approximate solution to Ax = b
max_iter– the number of iterations performed before the tolerance was reached
tol– the residual after the final iteration
Note
In this case isn't considered residual vector - Iterative template routine – CG CG follows the algorithm described on p. 15 in the SIAM Templates book.

Definition at line 26 of file cgbet.h.

◆ CG() [2/2]

template<class Matrix , class Vector , class Preconditioner , class Real >
int CG ( Matrix &  A,
Vector x,
const Vector b,
Preconditioner &  M,
Vector residual,
int64_t &  max_iter,
Real &  tol,
const int  FromCurrent 
)

CG solves the symmetric positive definite linear system $ A x = b $ using the Conjugate Gradient method.

Returns
The return value indicates convergence within max_iter (input) iterations (0), or no convergence within max_iter iterations (1).
Upon successful return, output arguments have the following values:
Parameters
A– matrix of the system
b– vector of the system
M– preconditioner matrix
x– approximate solution to Ax = b
max_iter– the number of iterations performed before the tolerance was reached
tol– the residual after the final iteration
residual– residual vector (return)
FromCurrent– for type of operation (MultAdd) Iterative template routine – CG
CG follows the algorithm described on p. 15 in the SIAM Templates book.

Definition at line 31 of file cg.h.

References Dot(), Norm(), TPZMatrix< TVar >::Print(), test::res, and TPZExtractVal::val().

Referenced by TPZSandlerExtended::ComputeCapCoVertexTangent(), TPZSandlerExtended::ComputeCapTangent(), TPZSandlerExtended::ComputeCapVertexTangent(), TPZSandlerExtended::ComputeFailureTangent(), TPZSandlerExtended::JacobianCoVertex(), TPZSandlerExtended::Jacobianf1(), TPZSandlerExtended::Jacobianf2(), TPZSandlerExtended::Jacobianf2CoVertex(), TPZSandlerExtended::JacobianVertex(), TPZSandlerExtended::Res1(), TPZSandlerExtended::Res2(), TPZSandlerExtended::Res2CoVertex(), and TPZMatrix< STATE >::SolveCG().

◆ CGS()

template<class Matrix , class Vector , class Preconditioner , class Real >
int CGS ( const Matrix &  A,
Vector x,
const Vector b,
const Preconditioner &  M,
int64_t &  max_iter,
Real &  tol 
)

CGS solves the unsymmetric linear system $ Ax = b $ using the Conjugate Gradient Squared method.

Returns
The return value indicates convergence within max_iter (input) iterations (0), or no convergence within max_iter iterations (1).
Upon successful return, output arguments have the following values:
Parameters
A– matrix of the system
b– vector of the system
M– preconditioner matrix
x– approximate solution to Ax = b
max_iter– the number of iterations performed before the tolerance was reached
tol– the residual after the final iteration Iterative template routine – CGS
CGS follows the algorithm described on p. 26 of the SIAM Templates book.

Definition at line 24 of file cgs.h.

◆ CHEBY()

template<class Matrix , class Vector , class Preconditioner , class Real , class Type >
int CHEBY ( const Matrix &  A,
Vector x,
const Vector b,
const Preconditioner &  M,
int64_t &  max_iter,
Real &  tol,
Type  eigmin,
Type  eigmax 
)

CHEBY solves the symmetric positive definite linear system $ Ax = b $ using the Preconditioned Chebyshev Method.

Returns
The return value indicates convergence within max_iter (input) iterations (0), or no convergence within max_iter iterations (1).
Upon successful return, output arguments have the following values:
Parameters
A– matrix of the system
b– vector of the system
M– preconditioner matrix
x– approximate solution to Ax = b
max_iter– the number of iterations performed before the tolerance was reached
tol– the residual after the final iteration
eigmin– minimun eigenvalue type
eigmax– maximun eigenvalue type Iterative template routine – CHEBY
CHEBY follows the algorithm described on p. 30 of the SIAM Templates book.

Definition at line 29 of file cheby.h.

References pzgeom::tol.

◆ GMRES()

template<class Operator , class Vector , class Preconditioner , class Matrix , class Real >
int GMRES ( Operator &  A,
Vector x,
const Vector b,
Preconditioner &  M,
Matrix &  H,
int &  m,
int64_t &  max_iter,
Real &  tol,
Vector residual,
const int  FromCurrent 
)

GMRES solves the unsymmetric linear system Ax = b using the Generalized Minimum Residual method.

Returns
The return value indicates convergence within max_iter (input) iterations (0), or no convergence within max_iter iterations (1).
Upon successful return, output arguments have the following values:
Parameters
AMatrix of the system
bVector of the system
MPreconditioner matrix
HAuxiliar matrix (?)
mSize (?)
xApproximate solution to $ Ax = b $
max_iterThe number of iterations performed before the tolerance was reached
tolThe residual after the final iteration
residualResidual vector (return)
FromCurrentFor type of operation (MultAdd) Iterative template routine – GMRES
GMRES follows the algorithm described on p. 20 of the SIAM Templates book.

Definition at line 71 of file gmres.h.

References ApplyPlaneRotation(), Dot(), fabs, GeneratePlaneRotation(), Norm(), test::res, Update(), and TPZExtractVal::val().

Referenced by TPZMatrix< STATE >::SolveGMRES().

◆ QMR()

template<class Matrix , class Vector , class Preconditioner1 , class Preconditioner2 , class Real >
int QMR ( const Matrix &  A,
Vector x,
const Vector b,
const Preconditioner1 &  M1,
const Preconditioner2 &  M2,
int64_t &  max_iter,
Real &  tol 
)

QMR solves the unsymmetric linear system $ Ax = b $ using the Quasi-Minimal Residual method.

Returns
Upon successful return, output arguments have the following values:
Parameters
AMatrix of the system
bVector of the system
M1Preconditioner matrix (first)
M2Preconditioner matrix (second)
xApproximate solution to $ Ax=b $
max_iterThe number of iterations performed before the tolerance was reached
tolThe residual after the final iteration (tolerance) Iterative template routine – QMR

Quasi-Minimal Residual method following the algorithm as described on p. 24 in the SIAM Templates book.

return value - indicates


0 - convergence within max_iter iterations 1 - no convergence after max_iter iterations breakdown in: 2 - rho 3 - beta 4 - gamma 5 - delta 6 - ep

7 - xi

Definition at line 40 of file qmr.h.

References gamma(), and sqrt.

◆ Update()

template<class Matrix , class Vector >
void Update ( Vector x,
int64_t  k,
Matrix &  h,
Vector s,
Vector  v[] 
)

Computes back solve. Updates on v vector.

Definition at line 22 of file gmres.h.

References h.

Referenced by GMRES().