Data Structures | Public Member Functions | Private Attributes | Friends
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...

Collaboration diagram for lduMatrix:
Collaboration graph
[legend]

Data Structures

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 &)
 Construct given an LDU addressed mesh and an Istream. More...
 
 ~lduMatrix ()
 Destructor. More...
 
const lduMeshmesh () const
 Return the LDU mesh from which 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 (scalarField &, const tmp< scalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
 Matrix multiplication with updated interfaces. More...
 
void Tmul (scalarField &, const tmp< scalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
 Matrix transpose multiplication with updated interfaces. More...
 
void sumA (scalarField &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &) const
 Sum the coefficients on each row of the matrix. More...
 
void residual (scalarField &rA, const scalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
 
tmp< scalarFieldresidual (const scalarField &psi, const scalarField &source, const FieldField< Field, scalar > &interfaceBouCoeffs, const lduInterfaceFieldPtrsList &interfaces, const direction cmpt) const
 
void initMatrixInterfaces (const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const scalarField &psiif, scalarField &result, const direction cmpt) const
 Initialise the update of interfaced interfaces. More...
 
void updateMatrixInterfaces (const FieldField< Field, scalar > &interfaceCoeffs, const lduInterfaceFieldPtrsList &interfaces, const scalarField &psiif, scalarField &result, const direction cmpt) const
 Update interfaced interfaces for matrix operations. 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
 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
 

Private Attributes

const lduMeshlduMesh_
 LDU mesh reference. More...
 
scalarFieldlowerPtr_
 Coefficients (not including interfaces) More...
 
scalarFielddiagPtr_
 
scalarFieldupperPtr_
 

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 hierachy starting from an empty matrix, then deriving diagonal, symmetric and asymmetric matrices.

Source files

Definition at line 77 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 40 of file lduMatrix.C.

◆ lduMatrix() [2/4]

lduMatrix ( const lduMatrix A)

Construct as copy.

Definition at line 49 of file lduMatrix.C.

References A(), lduMatrix::diagPtr_, lduMatrix::lowerPtr_, and lduMatrix::upperPtr_.

Here is the call graph for this function:

◆ lduMatrix() [3/4]

lduMatrix ( lduMatrix A,
bool  reUse 
)

Construct as copy or re-use as specified.

Definition at line 73 of file lduMatrix.C.

References A(), lduMatrix::diagPtr_, lduMatrix::lowerPtr_, and lduMatrix::upperPtr_.

Here is the call graph for this function:

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

References lduMatrix::diagPtr_, lduMatrix::hasDiag(), lduMatrix::lowerPtr_, and lduMatrix::upperPtr_.

Here is the call graph for this function:

◆ ~lduMatrix()

~lduMatrix ( )

Destructor.

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

References lduMatrix::lduMesh_.

Referenced by GAMGSolver::agglomerateMatrix(), algebraicPairGAMGAgglomeration::algebraicPairGAMGAgglomeration(), GAMGAgglomeration::New(), and GAMGSolver::procAgglomerateMatrix().

Here is the caller graph for this function:

◆ lduAddr()

const lduAddressing& lduAddr ( ) const
inline

◆ patchSchedule()

const lduSchedule& patchSchedule ( ) const
inline

Return the patch evaluation schedule.

Definition at line 551 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 212 of file lduMatrix.C.

◆ diag() [2/3]

Foam::scalarField & diag ( const label  nCoeffs)

Definition at line 230 of file lduMatrix.C.

◆ upper() [2/3]

Foam::scalarField & upper ( const label  nCoeffs)

Definition at line 241 of file lduMatrix.C.

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

References lduMatrix::diagPtr_.

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

Here is the caller graph for this function:

◆ hasUpper()

bool hasUpper ( ) const
inline

Definition at line 580 of file lduMatrix.H.

References lduMatrix::upperPtr_.

Referenced by Foam::operator<<(), and GAMGSolver::procAgglomerateMatrix().

Here is the caller graph for this function:

◆ hasLower()

bool hasLower ( ) const
inline

Definition at line 585 of file lduMatrix.H.

References lduMatrix::lowerPtr_.

Referenced by GAMGSolver::agglomerateMatrix(), Foam::operator<<(), and GAMGSolver::procAgglomerateMatrix().

Here is the caller graph for this function:

◆ diagonal()

bool diagonal ( ) const
inline

Definition at line 590 of file lduMatrix.H.

References lduMatrix::diagPtr_, lduMatrix::lowerPtr_, and lduMatrix::upperPtr_.

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

Here is the caller graph for this function:

◆ symmetric()

bool symmetric ( ) const
inline

Definition at line 595 of file lduMatrix.H.

References lduMatrix::diagPtr_, lduMatrix::lowerPtr_, and lduMatrix::upperPtr_.

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

References lduMatrix::diagPtr_, lduMatrix::lowerPtr_, and lduMatrix::upperPtr_.

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 50 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 68 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 ( scalarField Apsi,
const tmp< scalarField > &  tpsi,
const FieldField< Field, scalar > &  interfaceBouCoeffs,
const lduInterfaceFieldPtrsList interfaces,
const direction  cmpt 
) const

Matrix multiplication with updated interfaces.

Definition at line 35 of file lduMatrixATmul.C.

References tmp::clear(), Foam::diag(), and psi.

Here is the call graph for this function:

◆ Tmul()

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

Matrix transpose multiplication with updated interfaces.

Definition at line 96 of file lduMatrixATmul.C.

References tmp::clear(), Foam::diag(), and psi.

Here is the call graph for this function:

◆ sumA()

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

Sum the coefficients on each row of the matrix.

Definition at line 155 of file lduMatrixATmul.C.

References Foam::diag(), forAll, and UPtrList::set().

Here is the call graph for this function:

◆ residual() [1/2]

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

Definition at line 204 of file lduMatrixATmul.C.

References Foam::diag(), forAll, patchi, psi, and UPtrList::set().

Referenced by 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::scalarField > residual ( const scalarField psi,
const scalarField source,
const FieldField< Field, scalar > &  interfaceBouCoeffs,
const lduInterfaceFieldPtrsList interfaces,
const direction  cmpt 
) const

Definition at line 284 of file lduMatrixATmul.C.

References psi, and scalarField().

Here is the call graph for this function:

◆ initMatrixInterfaces()

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

Initialise the update of interfaced interfaces.

for matrix operations

Definition at line 31 of file lduMatrixUpdateMatrixInterfaces.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, UPtrList::set(), List::size(), and UPtrList::size().

Referenced by GAMGSolver::interpolate(), 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 FieldField< Field, scalar > &  interfaceCoeffs,
const lduInterfaceFieldPtrsList interfaces,
const scalarField psiif,
scalarField result,
const direction  cmpt 
) const

Update interfaced interfaces for matrix operations.

Definition at line 97 of file lduMatrixUpdateMatrixInterfaces.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, UPtrList::set(), List::size(), and UPtrList::size().

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

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

Return info proxy.

Used to print matrix information to a stream

Definition at line 704 of file lduMatrix.H.

◆ operator=()

void operator= ( const lduMatrix A)

Definition at line 88 of file lduMatrixOperations.C.

References A(), Foam::abort(), Foam::diag(), and Foam::FatalError.

Here is the call graph for this function:

◆ negate()

void negate ( )

Definition at line 125 of file lduMatrixOperations.C.

◆ operator+=()

void operator+= ( const lduMatrix A)

Definition at line 144 of file lduMatrixOperations.C.

References A(), Foam::diag(), Foam::endl(), Foam::nl, and WarningInFunction.

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

Here is the call graph for this function:

◆ operator*=() [1/2]

void operator*= ( const scalarField sf)

Definition at line 302 of file lduMatrixOperations.C.

References sf().

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

◆ H() [4/4]

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

Definition at line 69 of file lduMatrixTemplates.C.

◆ faceH() [3/4]

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

Definition at line 79 of file lduMatrixTemplates.C.

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

Field Documentation

◆ lduMesh_

const lduMesh& lduMesh_
private

LDU mesh reference.

Definition at line 82 of file lduMatrix.H.

Referenced by lduMatrix::lduAddr(), and lduMatrix::mesh().

◆ lowerPtr_

scalarField* lowerPtr_
private

Coefficients (not including interfaces)

Definition at line 85 of file lduMatrix.H.

Referenced by lduMatrix::asymmetric(), lduMatrix::diagonal(), lduMatrix::H(), lduMatrix::H1(), lduMatrix::hasLower(), lduMatrix::lduMatrix(), and lduMatrix::symmetric().

◆ diagPtr_

scalarField * diagPtr_
private

◆ upperPtr_

scalarField * upperPtr_
private

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