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...

Inheritance diagram for lduMatrix:
Collaboration diagram for lduMatrix:

Classes

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

Public Types

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

Public Member Functions

 ClassName ("lduMatrix")
 
 lduMatrix (const lduMesh &)
 Construct given an LDU addressed mesh. More...
 
 lduMatrix (const lduMatrix &)
 Construct as copy. More...
 
 lduMatrix (lduMatrix &, bool reuse)
 Construct as copy or re-use as specified. More...
 
 lduMatrix (const lduMesh &, Istream &)
 Construct given an LDU addressed mesh and an Istream from which the coefficients are read. More...
 
 ~lduMatrix ()
 Destructor. More...
 
const lduMeshmesh () const noexcept
 Return the LDU mesh from which the addressing is obtained. More...
 
void setLduMesh (const lduMesh &m)
 Set the LDU mesh containing the addressing is obtained. More...
 
const lduAddressinglduAddr () const
 Return the LDU addressing. More...
 
const lduSchedulepatchSchedule () const
 Return the patch evaluation schedule. More...
 
scalarFieldlower ()
 
scalarFielddiag ()
 
scalarFieldupper ()
 
scalarFieldlower (const label size)
 
scalarFielddiag (const label nCoeffs)
 
scalarFieldupper (const label nCoeffs)
 
const scalarFieldlower () const
 
const scalarFielddiag () const
 
const scalarFieldupper () const
 
bool hasDiag () const noexcept
 
bool hasUpper () const noexcept
 
bool hasLower () const noexcept
 
bool diagonal () const noexcept
 
bool symmetric () const noexcept
 
bool asymmetric () const noexcept
 
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. More...
 
void Tmul (solveScalarField &, const tmp< solveScalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
 Matrix transpose multiplication with updated interfaces. More...
 
void sumA (solveScalarField &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &) const
 Sum the coefficients on each row of the matrix. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
void operator= (const lduMatrix &)
 
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. More...
 
static constexpr const label defaultMaxIter = 1000
 Default maximum number of iterations for solvers (1000) More...
 
static const scalar defaultTolerance = 1e-6
 Default (absolute) tolerance (1e-6) More...
 

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 79 of file lduMatrix.H.

Member Enumeration Documentation

◆ normTypes

enum 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 103 of file lduMatrix.H.

Constructor & Destructor Documentation

◆ lduMatrix() [1/4]

lduMatrix ( const lduMesh mesh)

Construct given an LDU addressed mesh.

The coefficients are initially empty for subsequent setting.

Definition at line 54 of file lduMatrix.C.

◆ lduMatrix() [2/4]

lduMatrix ( const lduMatrix A)

Construct as copy.

Definition at line 63 of file lduMatrix.C.

References A.

◆ lduMatrix() [3/4]

lduMatrix ( lduMatrix A,
bool  reuse 
)

Construct as copy or re-use as specified.

Definition at line 87 of file lduMatrix.C.

References A.

◆ lduMatrix() [4/4]

lduMatrix ( const lduMesh mesh,
Istream is 
)

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

Definition at line 134 of file lduMatrix.C.

References lduMatrix::hasDiag().

Here is the call graph for this function:

◆ ~lduMatrix()

~lduMatrix ( )

Destructor.

Definition at line 160 of file lduMatrix.C.

Member Function Documentation

◆ ClassName()

ClassName ( "lduMatrix"  )

◆ mesh()

const lduMesh& mesh ( ) const
inlinenoexcept

Return the LDU mesh from which the addressing is obtained.

Definition at line 718 of file lduMatrix.H.

Referenced by assemblyFaceAreaPairGAMGAgglomeration::assemblyFaceAreaPairGAMGAgglomeration(), distributedDILUPreconditioner::distributedDILUPreconditioner(), lduMatrix::lduAddr(), and GAMGAgglomeration::New().

Here is the caller graph for this function:

◆ setLduMesh()

void setLduMesh ( const lduMesh m)
inline

Set the LDU mesh containing the addressing is obtained.

Definition at line 726 of file lduMatrix.H.

◆ lduAddr()

◆ patchSchedule()

const lduSchedule& patchSchedule ( ) const
inline

Return the patch evaluation schedule.

Definition at line 742 of file lduMatrix.H.

References lduMatrix::lduAddr(), and lduAddressing::patchSchedule().

Here is the call graph for this function:

◆ lower() [1/3]

◆ diag() [1/3]

◆ upper() [1/3]

◆ lower() [2/3]

Foam::scalarField & lower ( const label  size)

Definition at line 226 of file lduMatrix.C.

References Foam::Zero.

◆ diag() [2/3]

Foam::scalarField & diag ( const label  nCoeffs)

Definition at line 244 of file lduMatrix.C.

References Foam::Zero.

◆ upper() [2/3]

Foam::scalarField & upper ( const label  nCoeffs)

Definition at line 255 of file lduMatrix.C.

References Foam::Zero.

◆ lower() [3/3]

const Foam::scalarField & lower ( ) const

Definition at line 273 of file lduMatrix.C.

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

Here is the call graph for this function:

◆ diag() [3/3]

const Foam::scalarField & diag ( ) const

Definition at line 293 of file lduMatrix.C.

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

Here is the call graph for this function:

◆ upper() [3/3]

const Foam::scalarField & upper ( ) const

Definition at line 306 of file lduMatrix.C.

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

Here is the call graph for this function:

◆ hasDiag()

bool hasDiag ( ) const
inlinenoexcept

Definition at line 766 of file lduMatrix.H.

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

Here is the caller graph for this function:

◆ hasUpper()

bool hasUpper ( ) const
inlinenoexcept

Definition at line 771 of file lduMatrix.H.

Referenced by Foam::operator<<().

Here is the caller graph for this function:

◆ hasLower()

bool hasLower ( ) const
inlinenoexcept

Definition at line 776 of file lduMatrix.H.

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

Here is the caller graph for this function:

◆ diagonal()

bool diagonal ( ) const
inlinenoexcept

Definition at line 781 of file lduMatrix.H.

Referenced by lduMatrix::solver::New().

Here is the caller graph for this function:

◆ symmetric()

bool symmetric ( ) const
inlinenoexcept

Definition at line 786 of file lduMatrix.H.

Referenced by lduMatrix::solver::New(), lduMatrix::smoother::New(), and lduMatrix::preconditioner::New().

Here is the caller graph for this function:

◆ asymmetric()

◆ sumDiag()

void sumDiag ( )

◆ negSumDiag()

void negSumDiag ( )

◆ sumMagOffDiag()

void sumMagOffDiag ( scalarField sumOff) const

Definition at line 65 of file lduMatrixOperations.C.

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

Here is the call 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 32 of file lduMatrixATmul.C.

References UList< T >::begin(), tmp< T >::clear(), Foam::diag(), Foam::stringOps::lower(), psi, and Foam::stringOps::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 98 of file lduMatrixATmul.C.

References UList< T >::begin(), tmp< T >::clear(), Foam::diag(), Foam::stringOps::lower(), psi, and Foam::stringOps::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 162 of file lduMatrixATmul.C.

References UList< T >::begin(), Foam::diag(), forAll, Foam::stringOps::lower(), UPtrList< T >::set(), and Foam::stringOps::upper().

Here is the call 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

Definition at line 211 of file lduMatrixATmul.C.

References UList< T >::begin(), Foam::diag(), Foam::stringOps::lower(), psi, and Foam::stringOps::upper().

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

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

◆ 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 286 of file lduMatrixATmul.C.

References Foam::New(), and psi.

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 27 of file lduMatrixUpdateMatrixInterfaces.C.

References Foam::add(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, mesh, UPtrList< T >::set(), UList< T >::size(), and UPtrList< T >::size().

Referenced by GaussSeidelSmoother::smooth(), symGaussSeidelSmoother::smooth(), and nonBlockingGaussSeidelSmoother::smooth().

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

Update interfaced interfaces for matrix operations.

Definition at line 102 of file lduMatrixUpdateMatrixInterfaces.C.

References Foam::add(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, UPtrList< T >::get(), mesh, UPtrList< T >::set(), UPtrList< T >::size(), and lduInterfaceField::updateInterfaceMatrix().

Referenced by GaussSeidelSmoother::smooth(), symGaussSeidelSmoother::smooth(), and nonBlockingGaussSeidelSmoother::smooth().

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

◆ 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 327 of file lduMatrix.C.

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

Referenced by PCG::scalarSolve(), PBiCGStab::scalarSolve(), FPCG::scalarSolve(), PBiCG::solve(), and smoothSolver::solve().

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

◆ H() [1/4]

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

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

Here is the caller graph for this function:

◆ H() [2/4]

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

◆ H1()

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

Definition at line 300 of file lduMatrixATmul.C.

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

Referenced by fvMatrix< Type >::H1().

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

◆ faceH() [1/4]

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

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

Here is the caller graph for this function:

◆ faceH() [2/4]

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 920 of file lduMatrix.H.

◆ operator=()

void operator= ( const lduMatrix A)

Definition at line 85 of file lduMatrixOperations.C.

References A, Foam::diag(), Foam::stringOps::lower(), and Foam::stringOps::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:

◆ negate()

void negate ( )

Definition at line 119 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 138 of file lduMatrixOperations.C.

References A, Foam::ensightOutput::debug, Foam::diag(), Foam::endl(), Foam::stringOps::lower(), Foam::nl, Foam::stringOps::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 217 of file lduMatrixOperations.C.

References A, Foam::ensightOutput::debug, Foam::diag(), Foam::endl(), Foam::stringOps::lower(), Foam::nl, Foam::stringOps::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 296 of file lduMatrixOperations.C.

References Foam::stringOps::lower(), and Foam::stringOps::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*= ( scalar  s)

Definition at line 326 of file lduMatrixOperations.C.

References s.

◆ H() [3/4]

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

◆ H() [4/4]

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]

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, lduMatrix::lower(), Foam::New(), psi, UList< T >::size(), and lduMatrix::upper().

Here is the call graph for this function:

◆ faceH() [4/4]

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

Definition at line 104 of file lduMatrixTemplates.C.

Friends And Related Function Documentation

◆ operator<< [1/2]

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

◆ operator<< [2/2]

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

Member Data Documentation

◆ normTypesNames_

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

Names for the normTypes.

Definition at line 113 of file lduMatrix.H.

Referenced by LduMatrix< Type, DType, LUType >::solver::readControls(), and lduMatrix::solver::readControls().

◆ defaultMaxIter

constexpr const label defaultMaxIter = 1000
static

Default maximum number of iterations for solvers (1000)

Definition at line 118 of file lduMatrix.H.

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

◆ defaultTolerance

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

Default (absolute) tolerance (1e-6)

Definition at line 123 of file lduMatrix.H.

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


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