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 |
| 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 Types | |
| enum | normTypes : char { NO_NORM, DEFAULT_NORM, L1_SCALED_NORM } |
| Enumerated matrix normalisation types. More... | |
Public Member Functions | |
| ClassName ("lduMatrix") | |
| lduMatrix (const lduMesh &mesh) | |
| Construct (without coefficients) for an LDU addressed mesh. More... | |
| lduMatrix (const lduMatrix &) | |
| Copy construct. More... | |
| lduMatrix (lduMatrix &&) | |
| Move construct. More... | |
| lduMatrix (lduMatrix &, bool reuse) | |
| Construct as copy or re-use as specified. More... | |
| lduMatrix (const lduMesh &mesh, Istream &is) | |
| Construct given an LDU addressed mesh and an Istream from which the coefficients are read. More... | |
| ~lduMatrix ()=default | |
| Destructor. More... | |
| const lduMesh & | mesh () const noexcept |
| Return the LDU mesh from which the addressing is obtained. More... | |
| void | setLduMesh (const lduMesh &m) |
| Set the LDU mesh containing the addressing. More... | |
| const lduAddressing & | lduAddr () const |
| Return the LDU addressing. More... | |
| const lduSchedule & | patchSchedule () const |
| Return the patch evaluation schedule. More... | |
| const scalarField & | diag () const |
| const scalarField & | upper () const |
| const scalarField & | lower () const |
| scalarField & | diag () |
| scalarField & | upper () |
| scalarField & | lower () |
| scalarField & | diag (label size) |
| scalarField & | upper (label nCoeffs) |
| scalarField & | lower (label nCoeffs) |
| word | matrixTypeName () const |
| The matrix type (empty, diagonal, symmetric, ...) More... | |
| bool | hasDiag () const noexcept |
| bool | hasUpper () const noexcept |
| bool | hasLower () const noexcept |
| bool | diagonal () const noexcept |
| Matrix has diagonal only. More... | |
| bool | symmetric () const noexcept |
| Matrix is symmetric. More... | |
| bool | asymmetric () const noexcept |
| Matrix is asymmetric (ie, full) More... | |
| 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< solveScalarField > | residual (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 |
| Initialise the update of interfaced interfaces for matrix operations. More... | |
| 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 |
| Set the residual field using an IOField on the object registry if it exists. 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< scalarField > | H1 () 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< lduMatrix > | info () const noexcept |
| Return info proxy, used to print matrix information to a stream. More... | |
| void | operator= (const lduMatrix &) |
| Copy assignment. More... | |
| void | operator= (lduMatrix &&) |
| Move assignment. More... | |
| 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 |
Static Public Attributes | |
| static const Enum< normTypes > | normTypesNames_ |
| Names for the normTypes. More... | |
| static constexpr const label | defaultMaxIter = 1000 |
| Default maximum number of iterations for solvers (1000) More... | |
| static const scalar | defaultTolerance = 1e-6 |
| Default (absolute) tolerance (1e-6) More... | |
Friends | |
| Ostream & | operator<< (Ostream &, const lduMatrix &) |
| Ostream & | operator<< (Ostream &, const InfoProxy< lduMatrix > &) |
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.
Definition at line 80 of file lduMatrix.H.
|
strong |
Enumerated matrix normalisation types.
| Enumerator | |
|---|---|
| NO_NORM | "none" norm (returns 1) |
| DEFAULT_NORM | "default" norm (== L1_scaled) |
| L1_SCALED_NORM | "L1_scaled" norm |
Definition at line 114 of file lduMatrix.H.
Construct (without coefficients) for an LDU addressed mesh.
Not yet 'explicit' (legacy code may rely on implicit construct)
Definition at line 54 of file lduMatrix.C.
Move construct.
Definition at line 81 of file lduMatrix.C.
Construct given an LDU addressed mesh and an Istream from which the coefficients are read.
Definition at line 121 of file lduMatrix.C.
|
default |
Destructor.
| ClassName | ( | "lduMatrix" | ) |
|
inlinenoexcept |
Return the LDU mesh from which the addressing is obtained.
Definition at line 739 of file lduMatrix.H.
Referenced by assemblyFaceAreaPairGAMGAgglomeration::assemblyFaceAreaPairGAMGAgglomeration(), distributedDILUPreconditioner::distributedDILUPreconditioner(), lduMatrix::lduAddr(), GAMGAgglomeration::New(), and lduMatrix::patchSchedule().

|
inline |
Set the LDU mesh containing the addressing.
Definition at line 747 of file lduMatrix.H.
|
inline |
Return the LDU addressing.
Definition at line 755 of file lduMatrix.H.
References lduMesh::lduAddr(), and lduMatrix::mesh().
Referenced by DICPreconditioner::calcReciprocalD(), DILUPreconditioner::calcReciprocalD(), lduMatrix::H(), lduMatrix::H1(), cyclicFvPatchField< vector >::manipulateMatrix(), cyclicACMIFvPatchField< Type >::manipulateMatrix(), cyclicAMIFvPatchField< scalar >::manipulateMatrix(), symGaussSeidelSmoother::smooth(), GaussSeidelSmoother::smooth(), nonBlockingGaussSeidelSmoother::smooth(), and lduMatrix::sumDiag().


|
inline |
Return the patch evaluation schedule.
Definition at line 763 of file lduMatrix.H.
References lduMesh::lduAddr(), lduMatrix::mesh(), and lduAddressing::patchSchedule().

| const Foam::scalarField & diag | ( | ) | const |
Definition at line 163 of file lduMatrix.C.
References Foam::abort(), Foam::FatalError, and FatalErrorInFunction.
Referenced by velocityDampingConstraint::addDamping(), multiphaseInterSystem::addInterfacePorosity(), oversetFvMeshBase::addInterpolation(), interRegionExplicitPorositySource::addSup(), solidificationMeltingSource::addSup(), DICPreconditioner::DICPreconditioner(), DILUPreconditioner::DILUPreconditioner(), FDICPreconditioner::FDICPreconditioner(), EulerD2dt2Scheme< Type >::fvmD2dt2(), EulerDdtScheme< Type >::fvmDdt(), backwardDdtScheme< Type >::fvmDdt(), CoEulerDdtScheme< Type >::fvmDdt(), SLTSDdtScheme< Type >::fvmDdt(), localEulerDdtScheme< Type >::fvmDdt(), CrankNicolsonDdtScheme< Type >::fvmDdt(), GAMGSolver::GAMGSolver(), waWallFunctionFvPatchScalarField::manipulateMatrix(), mixedEnergyFvPatchScalarField::manipulateMatrix(), cyclicFvPatchField< vector >::manipulateMatrix(), oversetFvPatchField< Type >::manipulateMatrix(), cyclicACMIFvPatchField< Type >::manipulateMatrix(), cyclicAMIFvPatchField< scalar >::manipulateMatrix(), nonBlockingGaussSeidelSmoother::nonBlockingGaussSeidelSmoother(), oversetFvMeshBase::normalisation(), Foam::operator<<(), GaussSeidelSmoother::smooth(), symGaussSeidelSmoother::smooth(), nonBlockingGaussSeidelSmoother::smooth(), oversetFvMeshBase::solveOverset(), lduMatrix::sumDiag(), and oversetFvMeshBase::write().


| const Foam::scalarField & upper | ( | ) | const |
Definition at line 203 of file lduMatrix.C.
References Foam::abort(), Foam::FatalError, and FatalErrorInFunction.
Referenced by oversetFvMeshBase::addInterpolation(), algebraicPairGAMGAgglomeration::algebraicPairGAMGAgglomeration(), DICPreconditioner::calcReciprocalD(), DILUPreconditioner::calcReciprocalD(), lduMatrix::faceH(), oversetFvPatchField< Type >::fringeFlux(), gaussConvectionScheme< Type >::fvmDiv(), gaussLaplacianScheme< Type, GType >::fvmLaplacianUncorrected(), relaxedNonOrthoGaussLaplacianScheme< Type, GType >::fvmLaplacianUncorrected(), lduMatrix::H(), lduMatrix::H1(), waWallFunctionFvPatchScalarField::manipulateMatrix(), mixedEnergyFvPatchScalarField::manipulateMatrix(), cyclicFvPatchField< vector >::manipulateMatrix(), oversetFvPatchField< Type >::manipulateMatrix(), cyclicACMIFvPatchField< Type >::manipulateMatrix(), cyclicAMIFvPatchField< scalar >::manipulateMatrix(), lduMatrix::negSumDiag(), Foam::operator<<(), symGaussSeidelSmoother::smooth(), GaussSeidelSmoother::smooth(), nonBlockingGaussSeidelSmoother::smooth(), oversetFvMeshBase::solveOverset(), oversetFvPatchField< Type >::storeFringeCoefficients(), lduMatrix::sumDiag(), lduMatrix::sumMagOffDiag(), and oversetFvMeshBase::write().


| const Foam::scalarField & lower | ( | ) | const |
Definition at line 268 of file lduMatrix.C.
References Foam::abort(), Foam::FatalError, and FatalErrorInFunction.
Referenced by oversetFvMeshBase::addInterpolation(), algebraicPairGAMGAgglomeration::algebraicPairGAMGAgglomeration(), DILUPreconditioner::calcReciprocalD(), lduMatrix::faceH(), oversetFvPatchField< Type >::fringeFlux(), gaussConvectionScheme< Type >::fvmDiv(), GAMGSolver::GAMGSolver(), lduMatrix::H(), lduMatrix::H1(), waWallFunctionFvPatchScalarField::manipulateMatrix(), mixedEnergyFvPatchScalarField::manipulateMatrix(), cyclicFvPatchField< vector >::manipulateMatrix(), oversetFvPatchField< Type >::manipulateMatrix(), cyclicACMIFvPatchField< Type >::manipulateMatrix(), cyclicAMIFvPatchField< scalar >::manipulateMatrix(), lduMatrix::negSumDiag(), Foam::operator<<(), GaussSeidelSmoother::smooth(), symGaussSeidelSmoother::smooth(), nonBlockingGaussSeidelSmoother::smooth(), oversetFvMeshBase::solveOverset(), oversetFvPatchField< Type >::storeFringeCoefficients(), lduMatrix::sumDiag(), lduMatrix::sumMagOffDiag(), and oversetFvMeshBase::write().


| Foam::scalarField & diag | ( | ) |
Definition at line 176 of file lduMatrix.C.
| Foam::scalarField & upper | ( | ) |
Definition at line 223 of file lduMatrix.C.
| Foam::scalarField & lower | ( | ) |
Definition at line 288 of file lduMatrix.C.
| Foam::scalarField & diag | ( | label | size | ) |
Definition at line 188 of file lduMatrix.C.
| Foam::scalarField & upper | ( | label | nCoeffs | ) |
Definition at line 246 of file lduMatrix.C.
| Foam::scalarField & lower | ( | label | nCoeffs | ) |
Definition at line 311 of file lduMatrix.C.
| Foam::word matrixTypeName | ( | ) | const |
The matrix type (empty, diagonal, symmetric, ...)
Definition at line 146 of file lduMatrix.C.
|
inlinenoexcept |
Definition at line 794 of file lduMatrix.H.
Referenced by Foam::operator<<().

|
inlinenoexcept |
Definition at line 795 of file lduMatrix.H.
Referenced by Foam::operator<<().

|
inlinenoexcept |
Definition at line 796 of file lduMatrix.H.
Referenced by algebraicPairGAMGAgglomeration::algebraicPairGAMGAgglomeration(), and Foam::operator<<().

|
inlinenoexcept |
Matrix has diagonal only.
Definition at line 801 of file lduMatrix.H.
Referenced by lduMatrix::solver::New().

|
inlinenoexcept |
Matrix is symmetric.
Definition at line 809 of file lduMatrix.H.
Referenced by lduMatrix::solver::New(), lduMatrix::smoother::New(), and lduMatrix::preconditioner::New().

|
inlinenoexcept |
Matrix is asymmetric (ie, full)
Definition at line 817 of file lduMatrix.H.
Referenced by mixedEnergyFvPatchScalarField::manipulateMatrix(), cyclicFvPatchField< vector >::manipulateMatrix(), cyclicACMIFvPatchField< Type >::manipulateMatrix(), cyclicAMIFvPatchField< scalar >::manipulateMatrix(), lduMatrix::solver::New(), lduMatrix::smoother::New(), and lduMatrix::preconditioner::New().

| void sumDiag | ( | ) |
Definition at line 30 of file lduMatrixOperations.C.
References lduMatrix::diag(), lduMatrix::lduAddr(), lduMatrix::lower(), lduAddressing::lowerAddr(), UList< T >::size(), lduMatrix::upper(), and lduAddressing::upperAddr().

| void negSumDiag | ( | ) |
Definition at line 47 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().


| void sumMagOffDiag | ( | scalarField & | sumOff | ) | const |
Definition at line 65 of file lduMatrixOperations.C.
References lduMatrix::lower(), Foam::mag(), UList< T >::size(), and lduMatrix::upper().

| 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 32 of file lduMatrixATmul.C.
References UList< T >::begin(), tmp< T >::clear(), Foam::diag(), Foam::stringOps::lower(), psi, and Foam::stringOps::upper().

| 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 98 of file lduMatrixATmul.C.
References UList< T >::begin(), tmp< T >::clear(), Foam::diag(), Foam::stringOps::lower(), psi, and Foam::stringOps::upper().

| 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 162 of file lduMatrixATmul.C.
References UList< T >::begin(), Foam::diag(), forAll, Foam::stringOps::lower(), UPtrList< T >::set(), and Foam::stringOps::upper().

| 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 UList< T >::begin(), Foam::diag(), Foam::stringOps::lower(), psi, and Foam::stringOps::upper().
Referenced by faMatrix< Type >::residual(), and fvMatrix< Type >::residual().


| 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 Foam::New(), and psi.

| 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 27 of file lduMatrixUpdateMatrixInterfaces.C.
References Foam::add(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, mesh, UPtrList< T >::set(), UList< T >::size(), and UPtrList< T >::size().
Referenced by GaussSeidelSmoother::smooth(), symGaussSeidelSmoother::smooth(), and nonBlockingGaussSeidelSmoother::smooth().


| 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 102 of file lduMatrixUpdateMatrixInterfaces.C.
References Foam::add(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, UPtrList< T >::get(), mesh, UPtrList< T >::set(), UPtrList< T >::size(), and lduInterfaceField::updateInterfaceMatrix().
Referenced by GaussSeidelSmoother::smooth(), symGaussSeidelSmoother::smooth(), and nonBlockingGaussSeidelSmoother::smooth().


| 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 335 of file lduMatrix.C.
References DebugInfo, Foam::endl(), objectRegistry::findObject(), objectRegistry::getObjectPtr(), mesh, IOobject::scopedName(), and fvMesh::thisDb().
Referenced by PCG::scalarSolve(), PBiCGStab::scalarSolve(), FPCG::scalarSolve(), PBiCG::solve(), and smoothSolver::solve().


Referenced by faMatrix< Type >::H(), and fvMatrix< Type >::H().

| Foam::tmp< Foam::scalarField > H1 | ( | ) | const |
Definition at line 300 of file lduMatrixATmul.C.
References UList< T >::begin(), lduMatrix::lduAddr(), lduMatrix::lower(), lduAddressing::lowerAddr(), tmp< T >::New(), UList< T >::size(), lduMatrix::upper(), lduAddressing::upperAddr(), and Foam::Zero.
Referenced by fvMatrix< Type >::H1().


Referenced by faMatrix< Type >::flux(), and fvMatrix< Type >::flux().

Return info proxy, used to print matrix information to a stream.
Definition at line 946 of file lduMatrix.H.
| void operator= | ( | const lduMatrix & | A | ) |
Copy assignment.
Definition at line 85 of file lduMatrixOperations.C.
References A, Foam::diag(), Foam::stringOps::lower(), and Foam::stringOps::upper().
Referenced by faMatrix< Type >::operator=(), and fvMatrix< Type >::operator=().


| void operator= | ( | lduMatrix && | A | ) |
| void negate | ( | ) |
Definition at line 130 of file lduMatrixOperations.C.
Referenced by faMatrix< Type >::negate(), and fvMatrix< Type >::negate().

| void operator+= | ( | const lduMatrix & | A | ) |
Definition at line 149 of file lduMatrixOperations.C.
References A, Foam::ensightOutput::debug, Foam::diag(), Foam::endl(), Foam::stringOps::lower(), Foam::nl, Foam::stringOps::upper(), and WarningInFunction.
Referenced by faMatrix< Type >::operator+=(), and fvMatrix< Type >::operator+=().


| void operator-= | ( | const lduMatrix & | A | ) |
Definition at line 221 of file lduMatrixOperations.C.
References A, Foam::ensightOutput::debug, Foam::diag(), Foam::endl(), Foam::stringOps::lower(), Foam::nl, Foam::stringOps::upper(), and WarningInFunction.
Referenced by faMatrix< Type >::operator-=(), and fvMatrix< Type >::operator-=().


| void operator*= | ( | const scalarField & | sf | ) |
Definition at line 293 of file lduMatrixOperations.C.
References Foam::stringOps::lower(), and Foam::stringOps::upper().
Referenced by faMatrix< Type >::operator*=(), and fvMatrix< Type >::operator*=().


| void operator*= | ( | scalar | s | ) |
Definition at line 323 of file lduMatrixOperations.C.
References s.
| Foam::tmp<Foam::Field<Type> > H | ( | const Field< Type > & | psi | ) | const |
Definition at line 30 of file lduMatrixTemplates.C.
References UList< T >::begin(), lduMatrix::lduAddr(), lduMatrix::lower(), lduAddressing::lowerAddr(), Foam::New(), psi, UList< T >::size(), lduMatrix::upper(), and lduAddressing::upperAddr().

| Foam::tmp<Foam::Field<Type> > H | ( | const tmp< Field< Type >> & | tpsi | ) | const |
Definition at line 60 of file lduMatrixTemplates.C.
References H().

| Foam::tmp<Foam::Field<Type> > faceH | ( | const Field< Type > & | psi | ) | const |
Definition at line 70 of file lduMatrixTemplates.C.
References Foam::exit(), Foam::FatalError, FatalErrorInFunction, lduMatrix::lower(), Foam::New(), psi, UList< T >::size(), and lduMatrix::upper().

| Foam::tmp<Foam::Field<Type> > faceH | ( | const tmp< Field< Type >> & | tpsi | ) | const |
Definition at line 104 of file lduMatrixTemplates.C.
|
static |
Names for the normTypes.
Definition at line 124 of file lduMatrix.H.
Referenced by LduMatrix< Type, DType, LUType >::solver::readControls(), and lduMatrix::solver::readControls().
|
static |
Default maximum number of iterations for solvers (1000)
Definition at line 129 of file lduMatrix.H.
Referenced by lduMatrix::solver::readControls().
|
static |
Default (absolute) tolerance (1e-6)
Definition at line 134 of file lduMatrix.H.
Referenced by lduMatrix::solver::readControls().