Loading...
Searching...
No Matches
lduMatrix Class Reference

lduMatrix is a general matrix class in which the coefficients are stored as three arrays, one for the upper triangle, one for the lower triangle and a third for the diagonal. More...

#include <lduMatrix.H>

Inheritance diagram for lduMatrix:
Collaboration diagram for lduMatrix:

Classes

class  solver
 Abstract base-class for lduMatrix solvers. More...
class  smoother
 Abstract base-class for lduMatrix smoothers. More...
class  preconditioner
 Abstract base-class for lduMatrix preconditioners. More...

Public Types

enum class  normTypes : char { NO_NORM , DEFAULT_NORM , L1_SCALED_NORM }
 Enumerated matrix normalisation types. More...

Public Member Functions

 ClassName ("lduMatrix")
 lduMatrix (const lduMesh &mesh)
 Construct (without coefficients) for an LDU addressed mesh.
 lduMatrix (const lduMatrix &)
 Copy construct.
 lduMatrix (lduMatrix &&)
 Move construct.
 lduMatrix (lduMatrix &, bool reuse)
 Construct as copy or re-use as specified.
 lduMatrix (const lduMesh &mesh, Istream &is)
 Construct given an LDU addressed mesh and an Istream from which the coefficients are read.
 ~lduMatrix ()=default
 Destructor.
const lduMeshmesh () const noexcept
 Return the LDU mesh from which the addressing is obtained.
void setLduMesh (const lduMesh &m)
 Set the LDU mesh containing the addressing.
const lduAddressinglduAddr () const
 Return the LDU addressing.
const lduSchedulepatchSchedule () const
 Return the patch evaluation schedule.
const scalarFielddiag () const
const scalarFieldupper () const
const scalarFieldlower () const
scalarFielddiag ()
scalarFieldupper ()
scalarFieldlower ()
scalarFielddiag (label size)
scalarFieldupper (label nCoeffs)
scalarFieldlower (label nCoeffs)
bool hasLowerCSR () const noexcept
const scalarFieldlowerCSR () const
scalarFieldlowerCSR ()
const solveScalarFieldwork () const
 Work array.
solveScalarFieldwork (label size) const
 Work array.
word matrixTypeName () const
 The matrix type (empty, diagonal, symmetric, ...).
bool hasDiag () const noexcept
bool hasUpper () const noexcept
bool hasLower () const noexcept
bool diagonal () const noexcept
 Matrix has diagonal only.
bool symmetric () const noexcept
 Matrix is symmetric.
bool asymmetric () const noexcept
 Matrix is asymmetric (ie, full).
void sumDiag ()
void negSumDiag ()
void sumMagOffDiag (scalarField &sumOff) const
void Amul (solveScalarField &, const tmp< solveScalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
 Matrix multiplication with updated interfaces.
void Tmul (solveScalarField &, const tmp< solveScalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
 Matrix transpose multiplication with updated interfaces.
void sumA (solveScalarField &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &) const
 Sum the coefficients on each row of the matrix.
void residual (solveScalarField &rA, const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
tmp< solveScalarFieldresidual (const solveScalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
void initMatrixInterfaces (const bool add, const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const solveScalarField &psiif, solveScalarField &result, const direction cmpt) const
 Initialise the update of interfaced interfaces for matrix operations.
void updateMatrixInterfaces (const bool add, const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const solveScalarField &psiif, solveScalarField &result, const direction cmpt, const label startRequest) const
 Update interfaced interfaces for matrix operations.
void setResidualField (const scalarField &residual, const word &fieldName, const bool initial) const
 Set the residual field using an IOField on the object registry if it exists.
template<class Type>
tmp< Field< Type > > H (const Field< Type > &) const
template<class Type>
tmp< Field< Type > > H (const tmp< Field< Type > > &) const
tmp< scalarFieldH1 () const
template<class Type>
tmp< Field< Type > > faceH (const Field< Type > &) const
template<class Type>
tmp< Field< Type > > faceH (const tmp< Field< Type > > &) const
InfoProxy< lduMatrixinfo () const noexcept
 Return info proxy, used to print matrix information to a stream.
void operator= (const lduMatrix &)
 Copy assignment.
void operator= (lduMatrix &&)
 Move assignment.
void negate ()
void operator+= (const lduMatrix &)
void operator-= (const lduMatrix &)
void operator*= (const scalarField &)
void operator*= (scalar)
template<class Type>
Foam::tmp< Foam::Field< Type > > H (const Field< Type > &psi) const
template<class Type>
Foam::tmp< Foam::Field< Type > > H (const tmp< Field< Type > > &tpsi) const
template<class Type>
Foam::tmp< Foam::Field< Type > > faceH (const Field< Type > &psi) const
template<class Type>
Foam::tmp< Foam::Field< Type > > faceH (const tmp< Field< Type > > &tpsi) const

Static Public Attributes

static const Enum< normTypesnormTypesNames_
 Names for the normTypes.
static constexpr const label defaultMaxIter = 1000
 Default maximum number of iterations for solvers (1000).
static const scalar defaultTolerance = 1e-6
 Default (absolute) tolerance (1e-6).

Friends

Ostreamoperator<< (Ostream &, const lduMatrix &)
Ostreamoperator<< (Ostream &, const InfoProxy< lduMatrix > &)

Detailed Description

lduMatrix is a general matrix class in which the coefficients are stored as three arrays, one for the upper triangle, one for the lower triangle and a third for the diagonal.

Addressing arrays must be supplied for the upper and lower triangles.

It might be better if this class were organised as a hierarchy starting from an empty matrix, then deriving diagonal, symmetric and asymmetric matrices.

Source files

Definition at line 80 of file lduMatrix.H.

Member Enumeration Documentation

◆ normTypes

enum class normTypes : char
strong

Enumerated matrix normalisation types.

Enumerator
NO_NORM 

"none" norm (returns 1)

DEFAULT_NORM 

"default" norm (== L1_scaled)

L1_SCALED_NORM 

"L1_scaled" norm

Definition at line 124 of file lduMatrix.H.

Constructor & Destructor Documentation

◆ lduMatrix() [1/5]

◆ lduMatrix() [2/5]

lduMatrix ( const lduMatrix & A)

Copy construct.

Definition at line 60 of file lduMatrix.C.

References A, if(), and lduMatrix().

Here is the call graph for this function:

◆ lduMatrix() [3/5]

lduMatrix ( lduMatrix && A)

Move construct.

Definition at line 88 of file lduMatrix.C.

References A, and lduMatrix().

Here is the call graph for this function:

◆ lduMatrix() [4/5]

lduMatrix ( lduMatrix & A,
bool reuse )

Construct as copy or re-use as specified.

Definition at line 99 of file lduMatrix.C.

References A, and lduMatrix().

Here is the call graph for this function:

◆ lduMatrix() [5/5]

lduMatrix ( const lduMesh & mesh,
Istream & is )

Construct given an LDU addressed mesh and an Istream from which the coefficients are read.

Definition at line 139 of file lduMatrix.C.

References lowerCSR(), and mesh().

Here is the call graph for this function:

◆ ~lduMatrix()

~lduMatrix ( )
default

Destructor.

Member Function Documentation

◆ ClassName()

ClassName ( "lduMatrix" )

References lduMatrix(), and mesh().

Here is the call graph for this function:

◆ mesh()

◆ setLduMesh()

void setLduMesh ( const lduMesh & m)
inline

Set the LDU mesh containing the addressing.

Definition at line 761 of file lduMatrix.H.

Referenced by fvMatrix< scalar >::solveSegregated().

Here is the caller graph for this function:

◆ lduAddr()

◆ patchSchedule()

const lduSchedule & patchSchedule ( ) const
inline

Return the patch evaluation schedule.

Definition at line 777 of file lduMatrix.H.

References lduMesh::lduAddr(), mesh(), and lduAddressing::patchSchedule().

Referenced by initMatrixInterfaces(), and updateMatrixInterfaces().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ diag() [1/3]

const Foam::scalarField & diag ( ) const

Definition at line 195 of file lduMatrix.C.

References Foam::abort(), Foam::FatalError, and FatalErrorInFunction.

Referenced by faMatrix< Type >::addBoundaryDiag(), faMatrix< scalar >::addBoundaryDiag(), fvMatrix< Type >::addBoundaryDiag(), fvMatrix< scalar >::addBoundaryDiag(), faMatrix< Type >::addCmptAvBoundaryDiag(), faMatrix< scalar >::addCmptAvBoundaryDiag(), fvMatrix< Type >::addCmptAvBoundaryDiag(), fvMatrix< scalar >::addCmptAvBoundaryDiag(), velocityDampingConstraint::addDamping(), oversetFvMeshBase::addInterpolation(), interRegionExplicitPorositySource::addSup(), interRegionExplicitPorositySource::addSup(), solidificationMeltingSource::addSup(), Amul(), cyclicACMIFvPatchField< Type >::manipulateMatrix(), cyclicAMIFvPatchField< Type >::manipulateMatrix(), mixedEnergyFvPatchScalarField::manipulateMatrix(), oversetFvPatchField< Type >::manipulateMatrix(), waWallFunctionFvPatchScalarField::manipulateMatrix(), oversetFvMeshBase::normalisation(), Foam::operator<<(), faMatrix< Type >::setComponentReference(), faMatrix< scalar >::setComponentReference(), fvMatrix< Type >::setComponentReference(), fvMatrix< scalar >::setComponentReference(), fvMatrix< Type >::setValuesFromList(), faMatrix< Type >::solve(), faMatrix< scalar >::solve(), fvMatrix< Type >::solveCoupled(), oversetFvMeshBase::solveOverset(), fvMatrix< scalar >::solver(), fvMatrix< Type >::solveSegregated(), fvMatrix< scalar >::solveSegregated(), Foam::Expression::Sp(), sumDiag(), Foam::Expression::SuSp(), cyclicFvPatchField< scalar >::write(), and oversetFvMeshBase::write().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ upper() [1/3]

◆ lower() [1/3]

◆ diag() [2/3]

Foam::scalarField & diag ( )

Definition at line 208 of file lduMatrix.C.

References lduAddr().

Here is the call graph for this function:

◆ upper() [2/3]

Foam::scalarField & upper ( )

Definition at line 255 of file lduMatrix.C.

References lduAddr().

Here is the call graph for this function:

◆ lower() [2/3]

Foam::scalarField & lower ( )

Definition at line 326 of file lduMatrix.C.

References lduAddr().

Here is the call graph for this function:

◆ diag() [3/3]

Foam::scalarField & diag ( label size)

Definition at line 220 of file lduMatrix.C.

◆ upper() [3/3]

Foam::scalarField & upper ( label nCoeffs)

Definition at line 281 of file lduMatrix.C.

◆ lower() [3/3]

Foam::scalarField & lower ( label nCoeffs)

Definition at line 351 of file lduMatrix.C.

◆ hasLowerCSR()

bool hasLowerCSR ( ) const
inlinenoexcept

Definition at line 803 of file lduMatrix.H.

References Foam::noexcept.

Referenced by Amul(), and Foam::operator<<().

Here is the caller graph for this function:

◆ lowerCSR() [1/2]

const Foam::scalarField & lowerCSR ( ) const

Definition at line 376 of file lduMatrix.C.

References Foam::abort(), Foam::FatalError, FatalErrorInFunction, and lduAddr().

Referenced by Amul(), lduMatrix(), and Foam::operator<<().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ lowerCSR() [2/2]

Foam::scalarField & lowerCSR ( )

Definition at line 404 of file lduMatrix.C.

References Foam::abort(), Foam::FatalError, FatalErrorInFunction, and lduAddr().

Here is the call graph for this function:

◆ work() [1/2]

const Foam::solveScalarField & work ( ) const

Work array.

Definition at line 443 of file lduMatrix.C.

References Foam::abort(), Foam::FatalError, and FatalErrorInFunction.

Here is the call graph for this function:

◆ work() [2/2]

Foam::solveScalarField & work ( label size) const

Work array.

Definition at line 432 of file lduMatrix.C.

◆ matrixTypeName()

Foam::word matrixTypeName ( ) const

The matrix type (empty, diagonal, symmetric, ...).

Definition at line 178 of file lduMatrix.C.

Referenced by operator+=(), and operator-=().

Here is the caller graph for this function:

◆ hasDiag()

bool hasDiag ( ) const
inlinenoexcept

Definition at line 827 of file lduMatrix.H.

References Foam::noexcept.

Referenced by Foam::operator<<().

Here is the caller graph for this function:

◆ hasUpper()

bool hasUpper ( ) const
inlinenoexcept

Definition at line 828 of file lduMatrix.H.

References Foam::noexcept.

Referenced by Foam::operator<<().

Here is the caller graph for this function:

◆ hasLower()

bool hasLower ( ) const
inlinenoexcept

Definition at line 829 of file lduMatrix.H.

References Foam::noexcept.

Referenced by algebraicPairGAMGAgglomeration::algebraicPairGAMGAgglomeration(), and Foam::operator<<().

Here is the caller graph for this function:

◆ diagonal()

bool diagonal ( ) const
inlinenoexcept

Matrix has diagonal only.

Definition at line 834 of file lduMatrix.H.

References Foam::noexcept.

Referenced by operator+=(), and operator-=().

Here is the caller graph for this function:

◆ symmetric()

bool symmetric ( ) const
inlinenoexcept

Matrix is symmetric.

Definition at line 842 of file lduMatrix.H.

References Foam::noexcept.

Referenced by lduMatrix::preconditioner::New(), operator*=(), operator+=(), operator-=(), faMatrix< Type >::setValuesFromList(), and fvMatrix< Type >::setValuesFromList().

Here is the caller graph for this function:

◆ asymmetric()

◆ sumDiag()

void sumDiag ( )

Definition at line 30 of file lduMatrixOperations.C.

References diag(), lduAddr(), lduMatrix(), lduAddressing::lowerAddr(), UList< T >::size(), and lduAddressing::upperAddr().

Here is the call graph for this function:

◆ negSumDiag()

void negSumDiag ( )

Definition at line 47 of file lduMatrixOperations.C.

References Foam::diag(), lduAddr(), lduMatrix(), and UList< T >::size().

Here is the call graph for this function:

◆ sumMagOffDiag()

void sumMagOffDiag ( scalarField & sumOff) const

Definition at line 64 of file lduMatrixOperations.C.

References lduAddr(), lduMatrix(), Foam::mag(), and UList< T >::size().

Referenced by faMatrix< Type >::relax(), and fvMatrix< Type >::relax().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Amul()

void Amul ( solveScalarField & Apsi,
const tmp< solveScalarField > & tpsi,
const FieldField< Field, scalar > & interfaceBouCoeffs,
const lduInterfaceFieldPtrsList & interfaces,
const direction cmpt ) const

Matrix multiplication with updated interfaces.

Definition at line 31 of file lduMatrixATmul.C.

References UList< T >::begin(), tmp< T >::clear(), diag(), Foam::endl(), hasLowerCSR(), initMatrixInterfaces(), lduAddr(), lower(), lowerCSR(), UPstream::nRequests(), PoutInFunction, psi, UList< T >::size(), updateMatrixInterfaces(), and upper().

Here is the call graph for this function:

◆ Tmul()

void Tmul ( solveScalarField & Tpsi,
const tmp< solveScalarField > & tpsi,
const FieldField< Field, scalar > & interfaceIntCoeffs,
const lduInterfaceFieldPtrsList & interfaces,
const direction cmpt ) const

Matrix transpose multiplication with updated interfaces.

Definition at line 149 of file lduMatrixATmul.C.

References UList< T >::begin(), tmp< T >::clear(), Foam::diag(), initMatrixInterfaces(), lduAddr(), lower(), UPstream::nRequests(), psi, updateMatrixInterfaces(), and upper().

Here is the call graph for this function:

◆ sumA()

void sumA ( solveScalarField & sumA,
const FieldField< Field, scalar > & interfaceBouCoeffs,
const lduInterfaceFieldPtrsList & interfaces ) const

Sum the coefficients on each row of the matrix.

Definition at line 213 of file lduMatrixATmul.C.

References Foam::diag(), forAll, lduAddr(), lower(), UPtrList< T >::set(), sumA(), and upper().

Referenced by sumA().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ residual() [1/2]

void residual ( solveScalarField & rA,
const solveScalarField & psi,
const scalarField & source,
const FieldField< Field, scalar > & interfaceBouCoeffs,
const lduInterfaceFieldPtrsList & interfaces,
const direction cmpt ) const

◆ residual() [2/2]

Foam::tmp< Foam::Field< Foam::solveScalar > > residual ( const solveScalarField & psi,
const scalarField & source,
const FieldField< Field, scalar > & interfaceBouCoeffs,
const lduInterfaceFieldPtrsList & interfaces,
const direction cmpt ) const

Definition at line 337 of file lduMatrixATmul.C.

References tmp< T >::New(), psi, and residual().

Here is the call graph for this function:

◆ initMatrixInterfaces()

void initMatrixInterfaces ( const bool add,
const FieldField< Field, scalar > & interfaceCoeffs,
const lduInterfaceFieldPtrsList & interfaces,
const solveScalarField & psiif,
solveScalarField & result,
const direction cmpt ) const

Initialise the update of interfaced interfaces for matrix operations.

Definition at line 26 of file lduMatrixUpdateMatrixInterfaces.C.

References Foam::add(), UPstream::buffered, UPstream::commsTypeNames, UPstream::defaultCommsType, Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, lduAddr(), mesh(), UPstream::nonBlocking, patchSchedule(), UPstream::scheduled, UPtrList< T >::set(), and UPtrList< T >::size().

Referenced by Amul(), residual(), faMatrix< Type >::solve(), fvMatrix< Type >::solveSegregated(), and Tmul().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ updateMatrixInterfaces()

void updateMatrixInterfaces ( const bool add,
const FieldField< Field, scalar > & interfaceCoeffs,
const lduInterfaceFieldPtrsList & interfaces,
const solveScalarField & psiif,
solveScalarField & result,
const direction cmpt,
const label startRequest ) const

◆ setResidualField()

void setResidualField ( const scalarField & residual,
const word & fieldName,
const bool initial ) const

Set the residual field using an IOField on the object registry if it exists.

Definition at line 456 of file lduMatrix.C.

References DebugInfo, Foam::endl(), objectRegistry::findObject(), objectRegistry::getObjectPtr(), mesh, residual(), IOobject::scopedName(), and fvMesh::thisDb().

Referenced by GAMGSolver::solve().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ H() [1/4]

template<class Type>
tmp< Field< Type > > H ( const Field< Type > & ) const

Referenced by faMatrix< Type >::H(), faMatrix< scalar >::H(), fvMatrix< Type >::H(), and fvMatrix< scalar >::H().

Here is the caller graph for this function:

◆ H() [2/4]

template<class Type>
tmp< Field< Type > > H ( const tmp< Field< Type > > & ) const

◆ H1()

Foam::tmp< Foam::scalarField > H1 ( ) const

Definition at line 352 of file lduMatrixATmul.C.

References lduAddr(), lower(), tmp< T >::New(), upper(), and Foam::Zero.

Referenced by fvMatrix< Type >::H1(), and fvMatrix< scalar >::H1().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ faceH() [1/4]

template<class Type>
tmp< Field< Type > > faceH ( const Field< Type > & ) const

Referenced by faceH(), faMatrix< Type >::flux(), and fvMatrix< Type >::flux().

Here is the caller graph for this function:

◆ faceH() [2/4]

template<class Type>
tmp< Field< Type > > faceH ( const tmp< Field< Type > > & ) const

◆ info()

InfoProxy< lduMatrix > info ( ) const
inlinenoexcept

Return info proxy, used to print matrix information to a stream.

Definition at line 979 of file lduMatrix.H.

References Foam::noexcept.

◆ operator=() [1/2]

void operator= ( const lduMatrix & A)

Copy assignment.

Definition at line 85 of file lduMatrixOperations.C.

References A, Foam::diag(), lduMatrix(), lower(), and upper().

Referenced by faMatrix< Type >::operator=(), and fvMatrix< Type >::operator=().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=() [2/2]

void operator= ( lduMatrix && A)

Move assignment.

Definition at line 117 of file lduMatrixOperations.C.

References A, and lduMatrix().

Here is the call graph for this function:

◆ negate()

void negate ( )

Definition at line 130 of file lduMatrixOperations.C.

Referenced by faMatrix< Type >::negate(), and fvMatrix< Type >::negate().

Here is the caller graph for this function:

◆ operator+=()

void operator+= ( const lduMatrix & A)

Definition at line 149 of file lduMatrixOperations.C.

References A, asymmetric(), Foam::diag(), diagonal(), Foam::endl(), lduMatrix(), lower(), matrixTypeName(), Foam::nl, symmetric(), upper(), and WarningInFunction.

Referenced by faMatrix< Type >::operator+=(), and fvMatrix< Type >::operator+=().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator-=()

void operator-= ( const lduMatrix & A)

Definition at line 221 of file lduMatrixOperations.C.

References A, asymmetric(), Foam::diag(), diagonal(), Foam::endl(), lduMatrix(), lower(), matrixTypeName(), Foam::nl, symmetric(), upper(), and WarningInFunction.

Referenced by faMatrix< Type >::operator-=(), and fvMatrix< Type >::operator-=().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator*=() [1/2]

void operator*= ( const scalarField & sf)

Definition at line 293 of file lduMatrixOperations.C.

References asymmetric(), lduAddr(), lower(), symmetric(), and upper().

Referenced by faMatrix< Type >::operator*=(), faMatrix< Type >::operator*=(), fvMatrix< Type >::operator*=(), and fvMatrix< Type >::operator*=().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator*=() [2/2]

void operator*= ( scalar s)

Definition at line 323 of file lduMatrixOperations.C.

References s().

Here is the call graph for this function:

◆ H() [3/4]

template<class Type>
Foam::tmp< Foam::Field< Type > > H ( const Field< Type > & psi) const

Definition at line 30 of file lduMatrixTemplates.C.

References UList< T >::begin(), lduAddr(), lower(), lduAddressing::lowerAddr(), Foam::New(), psi, UList< T >::size(), upper(), and lduAddressing::upperAddr().

Here is the call graph for this function:

◆ H() [4/4]

template<class Type>
Foam::tmp< Foam::Field< Type > > H ( const tmp< Field< Type > > & tpsi) const

Definition at line 60 of file lduMatrixTemplates.C.

References H().

Here is the call graph for this function:

◆ faceH() [3/4]

template<class Type>
Foam::tmp< Foam::Field< Type > > faceH ( const Field< Type > & psi) const

Definition at line 70 of file lduMatrixTemplates.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, lduAddr(), lduMatrix(), Foam::New(), psi, and UList< T >::size().

Here is the call graph for this function:

◆ faceH() [4/4]

template<class Type>
Foam::tmp< Foam::Field< Type > > faceH ( const tmp< Field< Type > > & tpsi) const

Definition at line 104 of file lduMatrixTemplates.C.

References faceH().

Here is the call graph for this function:

◆ operator<< [1/2]

Ostream & operator<< ( Ostream & ,
const lduMatrix &  )
friend

References lduMatrix().

◆ operator<< [2/2]

Ostream & operator<< ( Ostream & ,
const InfoProxy< lduMatrix > &  )
friend

Member Data Documentation

◆ normTypesNames_

const Foam::Enum< Foam::lduMatrix::normTypes > normTypesNames_
static

◆ defaultMaxIter

const label defaultMaxIter = 1000
staticconstexpr

Default maximum number of iterations for solvers (1000).

Definition at line 139 of file lduMatrix.H.

Referenced by lduMatrix::solver::readControls(), PBiCG::solve(), and lduMatrix::solver::solver().

◆ defaultTolerance

const Foam::scalar defaultTolerance = 1e-6
static

Default (absolute) tolerance (1e-6).

Definition at line 144 of file lduMatrix.H.

Referenced by lduMatrix::solver::readControls(), and lduMatrix::solver::solver().


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