SolverPerformance is the class returned by the LduMatrix solver containing performance statistics. More...
#include <SolverPerformance.H>

Public Member Functions | |
| ClassName ("SolverPerformance") | |
| SolverPerformance () | |
| SolverPerformance (const word &solverName, const word &fieldName, const Type &iRes=pTraits< Type >::zero, const Type &fRes=pTraits< Type >::zero, const labelType &nIter=pTraits< labelType >::zero, const bool converged=false, const bool singular=false) | |
| const word & | solverName () const noexcept |
| Return solver name. | |
| word & | solverName () noexcept |
| Return solver name. | |
| const word & | fieldName () const noexcept |
| Return field name. | |
| const Type & | initialResidual () const noexcept |
| Return initial residual. | |
| Type & | initialResidual () noexcept |
| Return initial residual. | |
| const Type & | finalResidual () const noexcept |
| Return final residual. | |
| Type & | finalResidual () noexcept |
| Return final residual. | |
| const labelType & | nIterations () const noexcept |
| Return number of iterations. | |
| labelType & | nIterations () noexcept |
| Return number of iterations. | |
| bool | converged () const noexcept |
| Has the solver converged? | |
| bool | singular () const |
| Is the matrix singular? | |
| bool | checkConvergence (const Type &tolerance, const Type &relTolerance, const int logLevel=0) |
| Check, store and return convergence. | |
| bool | checkSingularity (const Type &residual) |
| Singularity test. | |
| void | print (Ostream &os) const |
| Print summary of solver performance to the given stream. | |
| void | replace (const label cmpt, const SolverPerformance< typename pTraits< Type >::cmptType > &sp) |
| Replace component based on the minimal SolverPerformance. | |
| SolverPerformance< typename pTraits< Type >::cmptType > | max () |
| Return the summary maximum of SolverPerformance<Type>. | |
| bool | operator!= (const SolverPerformance< Type > &) const |
| Foam::SolverPerformance< Foam::scalar > | max () |
| SolverPerformance< scalar > | max () |
Static Public Attributes | |
| static const scalar | great_ |
| Large Type for the use in solvers. | |
| static const scalar | small_ |
| Small Type for the use in solvers. | |
| static const scalar | vsmall_ |
| Very small Type for the use in solvers. | |
Friends | |
| SolverPerformance< Type > | Foam::max (const SolverPerformance< Type > &, const SolverPerformance< Type > &) |
| Return the element-wise maximum of two SolverPerformance<Type>s. | |
| Istream & | operator>> (Istream &, SolverPerformance< Type > &) |
| Ostream & | operator<< (Ostream &, const SolverPerformance< Type > &) |
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition at line 77 of file SolverPerformance.H.
|
inline |
Definition at line 118 of file SolverPerformance.H.
Referenced by max(), and operator!=().

|
inline |
Definition at line 128 of file SolverPerformance.H.
| ClassName | ( | "SolverPerformance< Type >" | ) |
|
inlinenoexcept |
Return solver name.
Definition at line 154 of file SolverPerformance.H.
Referenced by operator!=(), faMatrix< Type >::solve(), fvMatrix< Type >::solveSegregated(), and solverInfo::updateSolverInfo().

|
inlinenoexcept |
Return solver name.
Definition at line 162 of file SolverPerformance.H.
|
inlinenoexcept |
Return field name.
Definition at line 171 of file SolverPerformance.H.
Referenced by operator!=(), and meshState::setSolverPerformance().

|
inlinenoexcept |
Return initial residual.
Definition at line 180 of file SolverPerformance.H.
Referenced by advectionDiffusion::correct(), radiativeIntensityRay::correct(), age::execute(), operator!=(), FPCG::scalarSolve(), PBiCGStab::scalarSolve(), PCG::scalarSolve(), PPCG::scalarSolveCG(), sensitivitySurface::smoothSensitivities(), adjointEikonalSolver::solve(), populationBalanceModel::solve(), elasticityMotionSolver::solve(), GAMGSolver::solve(), laplacianMotionSolver::solve(), PBiCCCG< Type, DType, LUType >::solve(), PBiCG::solve(), PBiCICG< Type, DType, LUType >::solve(), PCICG< Type, DType, LUType >::solve(), pLaplacianMotionSolver::solve(), SmoothSolver< Type, DType, LUType >::solve(), smoothSolver::solve(), Helmholtz::solveEqn(), shapeDesignVariables::solveMeshMovementEqn(), and solverInfo::updateSolverInfo().

|
inlinenoexcept |
Return initial residual.
Definition at line 188 of file SolverPerformance.H.
|
inlinenoexcept |
Return final residual.
Definition at line 197 of file SolverPerformance.H.
Referenced by operator!=(), FPCG::scalarSolve(), PBiCGStab::scalarSolve(), PCG::scalarSolve(), PPCG::scalarSolveCG(), GAMGSolver::solve(), PBiCCCG< Type, DType, LUType >::solve(), PBiCG::solve(), PBiCICG< Type, DType, LUType >::solve(), PCICG< Type, DType, LUType >::solve(), SmoothSolver< Type, DType, LUType >::solve(), smoothSolver::solve(), and solverInfo::updateSolverInfo().

|
inlinenoexcept |
Return final residual.
Definition at line 205 of file SolverPerformance.H.
|
inlinenoexcept |
Return number of iterations.
Definition at line 214 of file SolverPerformance.H.
Referenced by operator!=(), FPCG::scalarSolve(), PBiCGStab::scalarSolve(), PCG::scalarSolve(), PPCG::scalarSolveCG(), GAMGSolver::solve(), PBiCCCG< Type, DType, LUType >::solve(), PBiCG::solve(), PBiCICG< Type, DType, LUType >::solve(), PCICG< Type, DType, LUType >::solve(), SmoothSolver< Type, DType, LUType >::solve(), smoothSolver::solve(), and solverInfo::updateSolverInfo().

|
inlinenoexcept |
Return number of iterations.
Definition at line 222 of file SolverPerformance.H.
|
inlinenoexcept |
Has the solver converged?
Definition at line 231 of file SolverPerformance.H.
Referenced by operator!=(), and solverInfo::updateSolverInfo().

| bool singular | ( | ) | const |
Is the matrix singular?
Definition at line 44 of file SolverPerformance.C.
Referenced by checkSingularity(), max(), and operator!=().

| bool checkConvergence | ( | const Type & | tolerance, |
| const Type & | relTolerance, | ||
| const int | logLevel = 0 ) |
Check, store and return convergence.
Definition at line 56 of file SolverPerformance.C.
References Foam::cmptMultiply(), Foam::endl(), Foam::Info, and small_.
Referenced by FPCG::scalarSolve(), PBiCGStab::scalarSolve(), PCG::scalarSolve(), PPCG::scalarSolveCG(), GAMGSolver::solve(), PBiCCCG< Type, DType, LUType >::solve(), PBiCG::solve(), PBiCICG< Type, DType, LUType >::solve(), PCICG< Type, DType, LUType >::solve(), SmoothSolver< Type, DType, LUType >::solve(), and smoothSolver::solve().


| bool checkSingularity | ( | const Type & | residual | ) |
Singularity test.
Definition at line 28 of file SolverPerformance.C.
References Foam::component(), singular(), and vsmall_.
Referenced by FPCG::scalarSolve(), PBiCGStab::scalarSolve(), PCG::scalarSolve(), PBiCCCG< Type, DType, LUType >::solve(), PBiCG::solve(), PBiCICG< Type, DType, LUType >::solve(), and PCICG< Type, DType, LUType >::solve().


| void print | ( | Ostream & | os | ) | const |
Print summary of solver performance to the given stream.
Definition at line 85 of file SolverPerformance.C.
References Foam::component(), Foam::endl(), and os().
Referenced by viewFactor::calculate(), faMatrix< Type >::solve(), faMatrix< scalar >::solve(), fvMatrix< Type >::fvSolver::solve(), GAMGSolver::solve(), fvMatrix< Type >::solveCoupled(), fvMatrix< Type >::solveSegregated(), and fvMatrix< scalar >::solveSegregated().


| void replace | ( | const label | cmpt, |
| const SolverPerformance< typename pTraits< Type >::cmptType > & | sp ) |
Replace component based on the minimal SolverPerformance.
Definition at line 118 of file SolverPerformance.C.
Referenced by faMatrix< Type >::solve(), and fvMatrix< Type >::solveSegregated().

| Foam::SolverPerformance< typename Foam::pTraits< Type >::cmptType > max | ( | ) |
Return the summary maximum of SolverPerformance<Type>.
Effectively it will mostly return solverPerformanceScalar
Definition at line 133 of file SolverPerformance.C.
References Foam::cmptMax(), singular(), and SolverPerformance().

| bool operator!= | ( | const SolverPerformance< Type > & | sp | ) | const |
Definition at line 149 of file SolverPerformance.C.
References converged(), fieldName(), finalResidual(), initialResidual(), nIterations(), singular(), solverName(), and SolverPerformance().

| Foam::SolverPerformance< Foam::scalar > max | ( | ) |
Definition at line 36 of file solverPerformance.C.
| SolverPerformance< scalar > max | ( | ) |
|
friend |
Return the element-wise maximum of two SolverPerformance<Type>s.
|
friend |
|
friend |
|
static |
Large Type for the use in solvers.
Definition at line 103 of file SolverPerformance.H.
Referenced by FPCG::scalarSolve(), PCG::scalarSolve(), PBiCICG< Type, DType, LUType >::solve(), and PCICG< Type, DType, LUType >::solve().
|
static |
Small Type for the use in solvers.
Definition at line 108 of file SolverPerformance.H.
Referenced by checkConvergence(), and LduMatrix< scalar, scalar, scalar >::solver::readControls().
|
static |
Very small Type for the use in solvers.
Definition at line 113 of file SolverPerformance.H.
Referenced by checkSingularity(), PBiCICG< Type, DType, LUType >::solve(), and PCICG< Type, DType, LUType >::solve().