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:
[legend]

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 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 &)
 
 ~lduMatrix ()
 Destructor. More...
 
const lduMeshmesh () const
 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
 
bool hasUpper () const
 
bool hasLower () const
 
bool diagonal () const
 
bool symmetric () const
 
bool asymmetric () const
 
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
 
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
 
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
 Return info proxy. 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
 

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 83 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 49 of file lduMatrix.C.

◆ lduMatrix() [2/4]

lduMatrix ( const lduMatrix A)

Construct as copy.

Definition at line 58 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 82 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 129 of file lduMatrix.C.

References lduMatrix::hasDiag().

Here is the call graph for this function:

◆ ~lduMatrix()

~lduMatrix ( )

Destructor.

Definition at line 155 of file lduMatrix.C.

Member Function Documentation

◆ ClassName()

ClassName ( "lduMatrix"  )

◆ mesh()

const lduMesh & mesh ( ) const
inline

Return the LDU mesh from which the addressing is obtained.

Definition at line 566 of file lduMatrix.H.

Referenced by algebraicPairGAMGAgglomeration::algebraicPairGAMGAgglomeration(), assemblyFaceAreaPairGAMGAgglomeration::assemblyFaceAreaPairGAMGAgglomeration(), lduMatrix::initMatrixInterfaces(), lduMatrix::lduAddr(), GAMGAgglomeration::New(), faMatrix< Type >::solve(), and fvMatrix< Type >::solveSegregated().

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

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

Here is the caller graph for this function:

◆ lduAddr()

◆ patchSchedule()

const lduSchedule & patchSchedule ( ) const
inline

Return the patch evaluation schedule.

Definition at line 584 of file lduMatrix.H.

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

Referenced by lduMatrix::initMatrixInterfaces().

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

◆ lower() [1/3]

◆ diag() [1/3]

Foam::scalarField & diag ( )

Definition at line 192 of file lduMatrix.C.

References Foam::Zero.

Referenced by velocityDampingConstraint::addDamping(), multiphaseInterSystem::addInterfacePorosity(), dynamicOversetFvMesh::addInterpolation(), solidificationMeltingSource::addSup(), interRegionExplicitPorositySource::addSup(), lduMatrix::Amul(), diagonalPreconditioner::diagonalPreconditioner(), DiagonalPreconditioner< Type, DType, LUType >::DiagonalPreconditioner(), DICPreconditioner::DICPreconditioner(), DICSmoother::DICSmoother(), DILUPreconditioner::DILUPreconditioner(), DILUSmoother::DILUSmoother(), FDICPreconditioner::FDICPreconditioner(), FDICSmoother::FDICSmoother(), EulerD2dt2Scheme< Type >::fvmD2dt2(), backwardDdtScheme< Type >::fvmDdt(), CoEulerDdtScheme< Type >::fvmDdt(), CrankNicolsonDdtScheme< Type >::fvmDdt(), EulerDdtScheme< Type >::fvmDdt(), localEulerDdtScheme< Type >::fvmDdt(), SLTSDdtScheme< Type >::fvmDdt(), GAMGSolver::GAMGSolver(), mixedEnergyFvPatchScalarField::manipulateMatrix(), waWallFunctionFvPatchScalarField::manipulateMatrix(), cyclicFvPatchField< Type >::manipulateMatrix(), cyclicACMIFvPatchField< Type >::manipulateMatrix(), cyclicAMIFvPatchField< Type >::manipulateMatrix(), nonBlockingGaussSeidelSmoother::nonBlockingGaussSeidelSmoother(), dynamicOversetFvMesh::normalisation(), Foam::operator<<(), nonBlockingGaussSeidelSmoother::smooth(), GaussSeidelSmoother::smooth(), symGaussSeidelSmoother::smooth(), fvMatrix< Type >::solveSegregated(), lduMatrix::sumDiag(), and dynamicOversetFvMesh::write().

Here is the caller graph for this function:

◆ upper() [1/3]

◆ lower() [2/3]

Foam::scalarField & lower ( const label  size)

Definition at line 221 of file lduMatrix.C.

References Foam::Zero.

◆ diag() [2/3]

Foam::scalarField & diag ( const label  nCoeffs)

Definition at line 239 of file lduMatrix.C.

References Foam::Zero.

◆ upper() [2/3]

Foam::scalarField & upper ( const label  nCoeffs)

Definition at line 250 of file lduMatrix.C.

References Foam::Zero.

◆ lower() [3/3]

const Foam::scalarField & lower ( ) const

Definition at line 268 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 288 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 301 of file lduMatrix.C.

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

Here is the call graph for this function:

◆ hasDiag()

bool hasDiag ( ) const
inline

Definition at line 608 of file lduMatrix.H.

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

Here is the caller graph for this function:

◆ hasUpper()

bool hasUpper ( ) const
inline

Definition at line 613 of file lduMatrix.H.

Referenced by Foam::operator<<().

Here is the caller graph for this function:

◆ hasLower()

bool hasLower ( ) const
inline

Definition at line 618 of file lduMatrix.H.

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

Here is the caller graph for this function:

◆ diagonal()

bool diagonal ( ) const
inline

Definition at line 623 of file lduMatrix.H.

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

Here is the caller graph for this function:

◆ symmetric()

bool symmetric ( ) const
inline

Definition at line 628 of file lduMatrix.H.

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

Here is the caller graph for this function:

◆ asymmetric()

bool asymmetric ( ) const
inline

◆ sumDiag()

void sumDiag ( )

Definition at line 36 of file lduMatrixOperations.C.

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

Here is the call graph for this function:

◆ negSumDiag()

void negSumDiag ( )

Definition at line 53 of file lduMatrixOperations.C.

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

Referenced by gaussConvectionScheme< Type >::fvmDiv(), gaussLaplacianScheme< Type, GType >::fvmLaplacianUncorrected(), and relaxedNonOrthoGaussLaplacianScheme< Type, GType >::fvmLaplacianUncorrected().

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

◆ sumMagOffDiag()

void sumMagOffDiag ( scalarField sumOff) const

Definition at line 70 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 37 of file lduMatrixATmul.C.

References UList< T >::begin(), tmp< T >::clear(), lduMatrix::diag(), lduMatrix::initMatrixInterfaces(), lduMatrix::lduAddr(), lduMatrix::lower(), lduAddressing::lowerAddr(), nFaces(), UPstream::nRequests(), psi, UList< T >::size(), lduMatrix::updateMatrixInterfaces(), lduMatrix::upper(), and lduAddressing::upperAddr().

Referenced by GAMGSolver::solve().

Here is the call graph for this function:
Here is the caller 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 103 of file lduMatrixATmul.C.

References UList< T >::begin(), tmp< T >::clear(), Foam::diag(), nFaces(), UPstream::nRequests(), and psi.

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

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

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

References UList< T >::begin(), Foam::diag(), nFaces(), UPstream::nRequests(), and psi.

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

References psi, tmp< T >::ref(), and UList< T >::size().

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

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

Referenced by lduMatrix::Amul(), nonBlockingGaussSeidelSmoother::smooth(), GaussSeidelSmoother::smooth(), and symGaussSeidelSmoother::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 106 of file lduMatrixUpdateMatrixInterfaces.C.

References Foam::add(), UPstream::blocking, UPstream::commsTypeNames, UPstream::defaultCommsType, Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, mesh, UPstream::nonBlocking, UPstream::nPollProcInterfaces, UPstream::parRun(), UPstream::resetRequests(), UPstream::scheduled, UPtrList< T >::set(), UPtrList< T >::size(), UList< T >::size(), and UPstream::waitRequests().

Referenced by lduMatrix::Amul(), nonBlockingGaussSeidelSmoother::smooth(), GaussSeidelSmoother::smooth(), and symGaussSeidelSmoother::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 321 of file lduMatrix.C.

References DebugInfo, Foam::endl(), objectRegistry::findObject(), dictionary::found(), objectRegistry::getObjectPtr(), mesh, 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]

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

◆ H() [2/4]

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

◆ H1()

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

Definition at line 306 of file lduMatrixATmul.C.

References Time::New(), nFaces(), 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
inline

Return info proxy.

Used to print matrix information to a stream

Definition at line 748 of file lduMatrix.H.

◆ operator=()

void operator= ( const lduMatrix A)

Definition at line 91 of file lduMatrixOperations.C.

References A, and Foam::diag().

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 125 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 144 of file lduMatrixOperations.C.

References A, Foam::diag(), Foam::endl(), Foam::nl, 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 223 of file lduMatrixOperations.C.

References A, Foam::diag(), Foam::endl(), Foam::nl, 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 302 of file lduMatrixOperations.C.

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

Here is the caller graph for this function:

◆ operator*=() [2/2]

void operator*= ( scalar  s)

Definition at line 332 of file lduMatrixOperations.C.

References s().

Here is the call graph for this function:

◆ H() [3/4]

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

Definition at line 36 of file lduMatrixTemplates.C.

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

Here is the call graph for this function:

◆ H() [4/4]

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

Definition at line 71 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 81 of file lduMatrixTemplates.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, lduMatrix::lower(), psi, tmp< T >::ref(), 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 115 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

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