DILUPreconditioner Class Reference

Simplified diagonal-based incomplete LU preconditioner for asymmetric matrices. The reciprocal of the preconditioned diagonal is calculated and stored. More...

Inheritance diagram for DILUPreconditioner:
[legend]
Collaboration diagram for DILUPreconditioner:
[legend]

Public Member Functions

 TypeName ("DILU")
 Runtime type information. More...
 
 DILUPreconditioner (const lduMatrix::solver &, const dictionary &solverControlsUnused)
 Construct from matrix components and preconditioner solver controls. More...
 
virtual ~DILUPreconditioner ()=default
 Destructor. More...
 
virtual void precondition (solveScalarField &wA, const solveScalarField &rA, const direction cmpt=0) const
 Return wA the preconditioned form of residual rA. More...
 
virtual void preconditionT (solveScalarField &wT, const solveScalarField &rT, const direction cmpt=0) const
 Return wT the transpose-matrix preconditioned form of residual rT. More...
 
- Public Member Functions inherited from lduMatrix::preconditioner
virtual const wordtype () const =0
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, preconditioner, symMatrix,(const solver &sol, const dictionary &solverControls),(sol, solverControls))
 
 declareRunTimeSelectionTable (autoPtr, preconditioner, asymMatrix,(const solver &sol, const dictionary &solverControls),(sol, solverControls))
 
 preconditioner (const solver &sol)
 
virtual ~preconditioner ()=default
 Destructor. More...
 
virtual void read (const dictionary &)
 

Static Public Member Functions

static void calcReciprocalD (solveScalarField &, const lduMatrix &)
 Calculate the reciprocal of the preconditioned diagonal. More...
 
- Static Public Member Functions inherited from lduMatrix::preconditioner
static word getName (const dictionary &)
 Find the preconditioner name (directly or from a sub-dictionary) More...
 
static autoPtr< preconditionerNew (const solver &sol, const dictionary &solverControls)
 Return a new preconditioner. More...
 

Additional Inherited Members

- Protected Attributes inherited from lduMatrix::preconditioner
const solversolver_
 Reference to the base-solver this preconditioner is used with. More...
 

Detailed Description

Simplified diagonal-based incomplete LU preconditioner for asymmetric matrices. The reciprocal of the preconditioned diagonal is calculated and stored.

Source files

Definition at line 56 of file DILUPreconditioner.H.

Constructor & Destructor Documentation

◆ DILUPreconditioner()

DILUPreconditioner ( const lduMatrix::solver sol,
const dictionary solverControlsUnused 
)

Construct from matrix components and preconditioner solver controls.

Definition at line 47 of file DILUPreconditioner.C.

References Foam::diag(), lduMatrix::diag(), and lduMatrix::solver::matrix().

Here is the call graph for this function:

◆ ~DILUPreconditioner()

virtual ~DILUPreconditioner ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName()

TypeName ( "DILU"  )

Runtime type information.

◆ calcReciprocalD()

void calcReciprocalD ( solveScalarField rD,
const lduMatrix matrix 
)
static

Calculate the reciprocal of the preconditioned diagonal.

Definition at line 65 of file DILUPreconditioner.C.

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

Here is the call graph for this function:

◆ precondition()

void precondition ( solveScalarField wA,
const solveScalarField rA,
const direction  cmpt = 0 
) const
virtual

Return wA the preconditioned form of residual rA.

Implements lduMatrix::preconditioner.

Definition at line 96 of file DILUPreconditioner.C.

◆ preconditionT()

void preconditionT ( solveScalarField wT,
const solveScalarField rT,
const direction  cmpt = 0 
) const
virtual

Return wT the transpose-matrix preconditioned form of residual rT.

Reimplemented from lduMatrix::preconditioner.

Definition at line 143 of file DILUPreconditioner.C.


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