OpenFOAM logo
The Open Source CFD Toolbox
  Source Guide OpenCFD Solutions Contact OpenFOAM

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:

List of all members.


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...
class  solverPerformance
 Class returned by the solver, containing performance statistics. More...

Public Member Functions

 ClassName ("lduMatrix")
 lduMatrix (const lduMesh &)
 Construct given an LDU addressed mesh.
 lduMatrix (const lduMatrix &)
 Construct as copy.
 lduMatrix (lduMatrix &, bool reUse)
 Construct as copy or re-use as specified.
 lduMatrix (const lduMesh &, Istream &)
 Construct given an LDU addressed mesh and an Istream.
 ~lduMatrix ()
const lduMeshmesh () const
 Return the LDU mesh from which the addressing is obtained.
const lduAddressinglduAddr () const
 Return the LDU addressing.
const lduSchedulepatchSchedule () const
 Return the patch evaluation schedule.
scalarFieldlower ()
scalarFielddiag ()
scalarFieldupper ()
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.
void Tmul (scalarField &, const tmp< scalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const
 Matrix transpose multiplication with updated interfaces.
void sumA (scalarField &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &) const
 Sum the coefficients on each row of the matrix.
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.
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.
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
void operator= (const lduMatrix &)
void negate ()
void operator+= (const lduMatrix &)
void operator-= (const lduMatrix &)
void operator*= (const scalarField &)
void operator*= (scalar)

Static Public Attributes

static const scalar great_ = 1.0e+20
 Large scalar for the use in solvers.
static const scalar small_ = 1.0e-20
 Small scalar for the use in solvers.

Friends

Ostreamoperator<< (Ostream &, const 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 71 of file lduMatrix.H.


Constructor & Destructor Documentation

lduMatrix ( const lduMesh mesh  ) 

Construct given an LDU addressed mesh.

The coefficients are initially empty for subsequent setting.

Definition at line 34 of file lduMatrix.C.

lduMatrix ( const lduMatrix A  ) 

Construct as copy.

Definition at line 43 of file lduMatrix.C.

lduMatrix ( lduMatrix A,
bool  reUse 
)

Construct as copy or re-use as specified.

Definition at line 67 of file lduMatrix.C.

lduMatrix ( const lduMesh mesh,
Istream is 
)

Construct given an LDU addressed mesh and an Istream.

from which the coefficients are read

Definition at line 115 of file lduMatrix.C.

~lduMatrix (  ) 

Definition at line 127 of file lduMatrix.C.


Member Function Documentation

ClassName ( "lduMatrix"   ) 

const lduMesh& mesh (  )  const [inline]

Return the LDU mesh from which the addressing is obtained.

Definition at line 676 of file lduMatrix.H.

References lduMesh::lduAddr().

Here is the call graph for this function:

const lduAddressing& lduAddr (  )  const [inline]

Return the LDU addressing.

Definition at line 682 of file lduMatrix.H.

References lduAddressing::patchSchedule().

Referenced by lduMatrix::Amul(), lduMatrix::diag(), lduMatrix::H(), lduMatrix::lower(), lduMatrix::sumDiag(), and lduMatrix::upper().

Here is the call graph for this function:

Here is the caller graph for this function:

const lduSchedule& patchSchedule (  )  const [inline]

Return the patch evaluation schedule.

Definition at line 688 of file lduMatrix.H.

Referenced by lduMatrix::initMatrixInterfaces().

Here is the caller graph for this function:

Foam::scalarField & lower (  ) 

Definition at line 146 of file lduMatrix.C.

Referenced by lduMatrix::Amul(), and lduMatrix::H().

Here is the caller graph for this function:

Foam::scalarField & diag (  ) 

Definition at line 164 of file lduMatrix.C.

References lduMatrix::lduAddr().

Referenced by lduMatrix::Amul(), fvMatrix< Type >::DD(), porousZone::modifyDdt(), lduMatrix::negate(), fvMatrix< Type >::relax(), fvMatrix< Type >::setComponentReference(), and lduMatrix::sumDiag().

Here is the call graph for this function:

Here is the caller graph for this function:

Foam::scalarField & upper (  ) 

Definition at line 175 of file lduMatrix.C.

References lduMatrix::lduAddr().

Referenced by lduMatrix::Amul(), and lduMatrix::H().

Here is the call graph for this function:

Here is the caller graph for this function:

const Foam::scalarField & lower (  )  const

Definition at line 193 of file lduMatrix.C.

References lduMatrix::lduAddr().

Here is the call graph for this function:

const Foam::scalarField & diag (  )  const

Definition at line 213 of file lduMatrix.C.

const Foam::scalarField & upper (  )  const

Definition at line 226 of file lduMatrix.C.

bool hasDiag (  )  const [inline]

Definition at line 704 of file lduMatrix.H.

bool hasUpper (  )  const [inline]

Definition at line 709 of file lduMatrix.H.

bool hasLower (  )  const [inline]

Definition at line 714 of file lduMatrix.H.

bool diagonal (  )  const [inline]

Definition at line 719 of file lduMatrix.H.

bool symmetric (  )  const [inline]

Definition at line 724 of file lduMatrix.H.

bool asymmetric (  )  const [inline]

Definition at line 729 of file lduMatrix.H.

void sumDiag (  ) 

Definition at line 26 of file lduMatrixOperations.C.

References lduMatrix::diag(), lduMatrix::lduAddr(), lduAddressing::lowerAddr(), and lduAddressing::upperAddr().

Here is the call graph for this function:

void negSumDiag (  ) 

Definition at line 43 of file lduMatrixOperations.C.

void sumMagOffDiag ( scalarField sumOff  )  const

Definition at line 61 of file lduMatrixOperations.C.

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

References UList< T >::begin(), lduMatrix::diag(), lduMatrix::initMatrixInterfaces(), lduMatrix::lduAddr(), lduMatrix::lower(), lduAddressing::lowerAddr(), psi, 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:

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

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

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

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

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

References Pstream::blocking, Pstream::defaultCommsType, forAll, Pstream::nonBlocking, lduMatrix::patchSchedule(), Pstream::scheduled, and List< T >::size().

Referenced by lduMatrix::Amul().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Referenced by lduMatrix::Amul().

Here is the caller graph for this function:

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

Definition at line 27 of file lduMatrixTemplates.C.

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

Here is the call graph for this function:

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

Definition at line 62 of file lduMatrixTemplates.C.

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

Reimplemented in fvMatrix< Type >, fvMatrix< Type >, and fvMatrix< Type >.

Definition at line 325 of file lduMatrixOperations.C.

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

Definition at line 72 of file lduMatrixTemplates.C.

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

Definition at line 106 of file lduMatrixTemplates.C.

void operator= ( const lduMatrix A  ) 

Definition at line 81 of file lduMatrixOperations.C.

void negate (  ) 

Reimplemented in fvMatrix< Type >.

Definition at line 118 of file lduMatrixOperations.C.

References lduMatrix::diag().

Here is the call graph for this function:

void operator+= ( const lduMatrix A  ) 

Definition at line 137 of file lduMatrixOperations.C.

void operator-= ( const lduMatrix A  ) 

Definition at line 205 of file lduMatrixOperations.C.

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

Here is the call graph for this function:

void operator*= ( const scalarField sf  ) 

Definition at line 273 of file lduMatrixOperations.C.

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

Here is the call graph for this function:

void operator*= ( scalar  s  ) 

Definition at line 306 of file lduMatrixOperations.C.


Friends And Related Function Documentation

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


Member Data Documentation

const Foam::scalar great_ = 1.0e+20 [static]

Large scalar for the use in solvers.

Definition at line 639 of file lduMatrix.H.

const Foam::scalar small_ = 1.0e-20 [static]

Small scalar for the use in solvers.

Definition at line 642 of file lduMatrix.H.


The documentation for this class was generated from the following files:
Copyright © 2000-2009 OpenCFD Ltd