Classes | Public Member Functions | Friends | List of all members
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...

Classes

class  preconditioner
 
class  smoother
 
class  solver
 

Public Member Functions

 ClassName ("lduMatrix")
 
 lduMatrix (const lduMesh &)
 
 lduMatrix (const lduMatrix &)
 
 lduMatrix (lduMatrix &, bool reuse)
 
 lduMatrix (const lduMesh &, Istream &)
 
 ~lduMatrix ()
 
const lduMeshmesh () const
 
void setLduMesh (const lduMesh &m)
 
const lduAddressinglduAddr () const
 
const lduSchedulepatchSchedule () const
 
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
 
void Tmul (solveScalarField &, const tmp< solveScalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
 
void sumA (solveScalarField &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &) const
 
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
 
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
 
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 79 of file lduMatrix.H.

Constructor & Destructor Documentation

◆ lduMatrix() [1/4]

lduMatrix ( const lduMesh mesh)

Definition at line 42 of file lduMatrix.C.

◆ lduMatrix() [2/4]

lduMatrix ( const lduMatrix A)

Definition at line 51 of file lduMatrix.C.

References A.

◆ lduMatrix() [3/4]

lduMatrix ( lduMatrix A,
bool  reuse 
)

Definition at line 75 of file lduMatrix.C.

References A.

◆ lduMatrix() [4/4]

lduMatrix ( const lduMesh mesh,
Istream is 
)

Definition at line 122 of file lduMatrix.C.

References lduMatrix::hasDiag().

Here is the call graph for this function:

◆ ~lduMatrix()

~lduMatrix ( )

Definition at line 148 of file lduMatrix.C.

Member Function Documentation

◆ ClassName()

ClassName ( "lduMatrix"  )

◆ mesh()

const lduMesh& mesh ( ) const
inline

◆ setLduMesh()

void setLduMesh ( const lduMesh m)
inline

Definition at line 568 of file lduMatrix.H.

◆ lduAddr()

const lduAddressing& lduAddr ( ) const
inline

◆ patchSchedule()

const lduSchedule& patchSchedule ( ) const
inline

Definition at line 580 of file lduMatrix.H.

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

Here is the call graph for this function:

◆ lower() [1/3]

const Foam::scalarField & lower ( )

◆ diag() [1/3]

const Foam::scalarField & diag ( )

◆ upper() [1/3]

const Foam::scalarField & upper ( )

◆ lower() [2/3]

Foam::scalarField & lower ( const label  size)

Definition at line 214 of file lduMatrix.C.

References Foam::Zero.

◆ diag() [2/3]

Foam::scalarField & diag ( const label  nCoeffs)

Definition at line 232 of file lduMatrix.C.

References Foam::Zero.

◆ upper() [2/3]

Foam::scalarField & upper ( const label  nCoeffs)

Definition at line 243 of file lduMatrix.C.

References Foam::Zero.

◆ lower() [3/3]

const scalarField& lower ( ) const

◆ diag() [3/3]

const scalarField& diag ( ) const

◆ upper() [3/3]

const scalarField& upper ( ) const

◆ hasDiag()

bool hasDiag ( ) const
inline

Definition at line 604 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 609 of file lduMatrix.H.

Referenced by Foam::operator<<().

Here is the caller graph for this function:

◆ hasLower()

bool hasLower ( ) const
inline

Definition at line 614 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 619 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 624 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()

bool asymmetric ( ) const
inline

Definition at line 629 of file lduMatrix.H.

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

Here is the caller graph for this function:

◆ sumDiag()

void sumDiag ( )

◆ negSumDiag()

void negSumDiag ( )

Definition at line 47 of file lduMatrixOperations.C.

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

Here is the call graph for this function:

◆ sumMagOffDiag()

void sumMagOffDiag ( scalarField sumOff) const

Definition at line 65 of file lduMatrixOperations.C.

References lduMatrix::lower(), Foam::mag(), UList::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

Definition at line 32 of file lduMatrixATmul.C.

References tmp::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

Definition at line 98 of file lduMatrixATmul.C.

References tmp::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

Definition at line 162 of file lduMatrixATmul.C.

References Foam::diag(), forAll, Foam::stringOps::lower(), UPtrList::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 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 psi, and tmp::ref().

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

Definition at line 27 of file lduMatrixUpdateMatrixInterfaces.C.

References Foam::add(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, mesh, UPtrList::set(), and UPtrList::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

Definition at line 100 of file lduMatrixUpdateMatrixInterfaces.C.

References Foam::add(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, mesh, UPtrList::set(), and UPtrList::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:

◆ setResidualField()

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

Definition at line 315 of file lduMatrix.C.

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

Referenced by PCG::scalarSolve(), PBiCGStab::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

◆ H() [2/4]

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

◆ H1()

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

◆ faceH() [1/4]

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

◆ faceH() [2/4]

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

◆ info()

InfoProxy<lduMatrix> info ( ) const
inline

Definition at line 744 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().

Here is the call graph for this function:

◆ negate()

void negate ( )

Definition at line 119 of file lduMatrixOperations.C.

◆ operator+=()

void operator+= ( const lduMatrix A)

◆ operator-=()

void operator-= ( const lduMatrix A)

◆ operator*=() [1/2]

void operator*= ( const scalarField sf)

Definition at line 296 of file lduMatrixOperations.C.

References Foam::stringOps::lower(), and Foam::stringOps::upper().

Here is the call 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 65 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 75 of file lduMatrixTemplates.C.

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