SolverPerformance is the class returned by the LduMatrix solver containing performance statistics. More...
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 |
Return solver name. More... | |
word & | solverName () |
Return solver name. More... | |
const word & | fieldName () const |
Return field name. More... | |
const Type & | initialResidual () const |
Return initial residual. More... | |
Type & | initialResidual () |
Return initial residual. More... | |
const Type & | finalResidual () const |
Return final residual. More... | |
Type & | finalResidual () |
Return final residual. More... | |
const labelType & | nIterations () const |
Return number of iterations. More... | |
labelType & | nIterations () |
Return number of iterations. More... | |
bool | converged () const |
Has the solver converged? More... | |
bool | singular () const |
Is the matrix singular? More... | |
bool | checkConvergence (const Type &tolerance, const Type &relTolerance) |
Check, store and return convergence. More... | |
bool | checkSingularity (const Type &residual) |
Singularity test. More... | |
void | print (Ostream &os) const |
Print summary of solver performance to the given stream. More... | |
void | replace (const label cmpt, const SolverPerformance< typename pTraits< Type >::cmptType > &sp) |
Replace component based on the minimal SolverPerformance. More... | |
SolverPerformance< typename pTraits< Type >::cmptType > | max () |
Return the summary maximum of SolverPerformance<Type> More... | |
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. More... | |
static const scalar | small_ |
Small Type for the use in solvers. More... | |
static const scalar | vsmall_ |
Very small Type for the use in solvers. More... | |
Friends | |
SolverPerformance< Type > | Foam::max (const SolverPerformance< Type > &, const SolverPerformance< Type > &) |
Return the element-wise maximum of two SolverPerformance<Type>s. More... | |
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 78 of file SolverPerformance.H.
|
inline |
Definition at line 113 of file SolverPerformance.H.
|
inline |
Definition at line 123 of file SolverPerformance.H.
ClassName | ( | "SolverPerformance< Type >" | ) |
|
inline |
Return solver name.
Definition at line 147 of file SolverPerformance.H.
Referenced by fvMatrix< Type >::solveSegregated().
|
inline |
Return solver name.
Definition at line 153 of file SolverPerformance.H.
|
inline |
Return field name.
Definition at line 160 of file SolverPerformance.H.
Referenced by Residuals< Type >::append().
|
inline |
Return initial residual.
Definition at line 167 of file SolverPerformance.H.
Referenced by radiativeIntensityRay::correct(), PBiCCCG< Type, DType, LUType >::solve(), PBiCICG< Type, DType, LUType >::solve(), PCICG< Type, DType, LUType >::solve(), SmoothSolver< Type, DType, LUType >::solve(), GAMGSolver::solve(), PBiCG::solve(), PBiCGStab::solve(), PCG::solve(), and smoothSolver::solve().
|
inline |
Return initial residual.
Definition at line 173 of file SolverPerformance.H.
|
inline |
Return final residual.
Definition at line 180 of file SolverPerformance.H.
Referenced by PBiCCCG< Type, DType, LUType >::solve(), PBiCICG< Type, DType, LUType >::solve(), PCICG< Type, DType, LUType >::solve(), SmoothSolver< Type, DType, LUType >::solve(), GAMGSolver::solve(), PBiCG::solve(), PBiCGStab::solve(), PCG::solve(), and smoothSolver::solve().
|
inline |
Return final residual.
Definition at line 186 of file SolverPerformance.H.
|
inline |
Return number of iterations.
Definition at line 193 of file SolverPerformance.H.
Referenced by PBiCCCG< Type, DType, LUType >::solve(), PBiCICG< Type, DType, LUType >::solve(), PCICG< Type, DType, LUType >::solve(), SmoothSolver< Type, DType, LUType >::solve(), GAMGSolver::solve(), PBiCG::solve(), PBiCGStab::solve(), PCG::solve(), and smoothSolver::solve().
|
inline |
Return number of iterations.
Definition at line 199 of file SolverPerformance.H.
|
inline |
Has the solver converged?
Definition at line 206 of file SolverPerformance.H.
bool singular |
Is the matrix singular?
Definition at line 48 of file SolverPerformance.C.
bool checkConvergence | ( | const Type & | tolerance, |
const Type & | relTolerance | ||
) |
Check, store and return convergence.
Definition at line 60 of file SolverPerformance.C.
References Foam::cmptMultiply(), Foam::endl(), and Foam::Info.
Referenced by PBiCCCG< Type, DType, LUType >::solve(), PBiCICG< Type, DType, LUType >::solve(), PCICG< Type, DType, LUType >::solve(), SmoothSolver< Type, DType, LUType >::solve(), GAMGSolver::solve(), PBiCG::solve(), PBiCGStab::solve(), PCG::solve(), and smoothSolver::solve().
bool checkSingularity | ( | const Type & | residual | ) |
Singularity test.
Definition at line 32 of file SolverPerformance.C.
References Foam::component().
Referenced by PBiCCCG< Type, DType, LUType >::solve(), PBiCICG< Type, DType, LUType >::solve(), PCICG< Type, DType, LUType >::solve(), PBiCG::solve(), PBiCGStab::solve(), and PCG::solve().
void print | ( | Ostream & | os | ) | const |
Print summary of solver performance to the given stream.
Definition at line 96 of file SolverPerformance.C.
References Foam::component(), Foam::endl(), and Foam::indent().
Referenced by fvMatrix< Type >::fvSolver::solve(), GAMGSolver::solve(), fvMatrix< Type >::solveCoupled(), and fvMatrix< Type >::solveSegregated().
void replace | ( | const label | cmpt, |
const SolverPerformance< typename pTraits< Type >::cmptType > & | sp | ||
) |
Replace component based on the minimal SolverPerformance.
Definition at line 130 of file SolverPerformance.C.
Foam::SolverPerformance< typename Foam::pTraits< Type >::cmptType > max |
Return the summary maximum of SolverPerformance<Type>
Effectively it will mostly return solverPerformanceScalar
Definition at line 145 of file SolverPerformance.C.
References Foam::cmptMax().
bool operator!= | ( | const SolverPerformance< Type > & | sp | ) | const |
Definition at line 161 of file SolverPerformance.C.
Foam::SolverPerformance< Foam::scalar > max | ( | ) |
Definition at line 41 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 102 of file SolverPerformance.H.
Referenced by PBiCICG< Type, DType, LUType >::solve(), PCICG< Type, DType, LUType >::solve(), and PCG::solve().
|
static |
Small Type for the use in solvers.
Definition at line 105 of file SolverPerformance.H.
Referenced by lduMatrix::solver::normFactor().
|
static |
Very small Type for the use in solvers.
Definition at line 108 of file SolverPerformance.H.
Referenced by PBiCICG< Type, DType, LUType >::solve(), and PCICG< Type, DType, LUType >::solve().