A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise. More...
#include <fvMatrix.H>


Classes | |
| class | fvSolver |
| Solver class returned by the solver function used for systems in which it is useful to cache the solver for reuse. More... | |
Public Types | |
| typedef GeometricField< Type, fvPatchField, volMesh > | psiFieldType |
| The geometric field type for psi. | |
| typedef GeometricField< Type, fvsPatchField, surfaceMesh > | faceFluxFieldType |
| Field type for face flux (for non-orthogonal correction). | |
| typedef std::unique_ptr< faceFluxFieldType > | faceFluxFieldPtrType |
| Declare return type of the faceFluxCorrectionPtr() function. | |
| Public Types inherited from lduMatrix | |
| enum class | normTypes : char { NO_NORM , DEFAULT_NORM , L1_SCALED_NORM } |
| Enumerated matrix normalisation types. More... | |
Public Member Functions | |
| ClassName ("fvMatrix") | |
| fvMatrix (const GeometricField< Type, fvPatchField, volMesh > &psi, const dimensionSet &ds) | |
| Construct given a field to solve for. | |
| fvMatrix (const fvMatrix< Type > &) | |
| Copy construct. | |
| fvMatrix (const tmp< fvMatrix< Type > > &) | |
| Copy/move construct from tmp<fvMatrix<Type>>. | |
| fvMatrix (const GeometricField< Type, fvPatchField, volMesh > &psi, Istream &is) | |
| Deprecated(2022-05) - construct with dimensionSet instead. | |
| tmp< fvMatrix< Type > > | clone () const |
| Construct and return a clone. | |
| virtual | ~fvMatrix () |
| Destructor. | |
| label | nMatrices () const |
| const fvMatrix< Type > & | matrix (const label i) const |
| fvMatrix< Type > & | matrix (const label i) |
| label | globalPatchID (const label fieldi, const label patchi) const |
| void | transferFvMatrixCoeffs () |
| Transfer lower, upper, diag and source to this fvMatrix. | |
| void | createOrUpdateLduPrimitiveAssembly () |
| Create or update ldu assembly. | |
| lduPrimitiveMeshAssembly * | lduMeshPtr () |
| Access to lduPrimitiveMeshAssembly. | |
| const lduPrimitiveMeshAssembly * | lduMeshPtr () const |
| Const Access to lduPrimitiveMeshAssembly. | |
| void | manipulateMatrix (direction cmp) |
| Manipulate matrix. | |
| void | setBounAndInterCoeffs () |
| Manipulate boundary/internal coeffs for coupling. | |
| void | setInterfaces (lduInterfaceFieldPtrsList &, PtrDynList< lduInterfaceField > &newInterfaces) |
| Set interfaces. | |
| void | mapContributions (label fieldi, const FieldField< Field, Type > &fluxContrib, FieldField< Field, Type > &contrib, bool internal) const |
| Add internal and boundary contribution to local patches. | |
| const lduPrimitiveMeshAssembly & | lduMeshAssembly () |
| Return optional lduAdressing. | |
| const GeometricField< Type, fvPatchField, volMesh > & | psi (const label i=0) const |
| Return psi. | |
| GeometricField< Type, fvPatchField, volMesh > & | psi (const label i=0) |
| void | clear () |
| Clear multiple fvMatrices. | |
| const dimensionSet & | dimensions () const noexcept |
| Field< Type > & | source () noexcept |
| const Field< Type > & | source () const noexcept |
| const FieldField< Field, Type > & | internalCoeffs () const noexcept |
| fvBoundary scalar field containing pseudo-matrix coeffs for internal cells | |
| FieldField< Field, Type > & | internalCoeffs () noexcept |
| fvBoundary scalar field containing pseudo-matrix coeffs for internal cells | |
| const FieldField< Field, Type > & | boundaryCoeffs () const noexcept |
| fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells | |
| FieldField< Field, Type > & | boundaryCoeffs () noexcept |
| fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells | |
| faceFluxFieldPtrType & | faceFluxCorrectionPtr () |
| Return pointer to face-flux non-orthogonal correction field. | |
| void | faceFluxCorrectionPtr (faceFluxFieldType *flux) |
| Set pointer to face-flux non-orthogonal correction field. | |
| bool | hasFaceFluxCorrection () const noexcept |
| True if face-flux non-orthogonal correction field exists. | |
| void | setValues (const labelUList &cellLabels, const Type &value) |
| Set solution in given cells to the specified value and eliminate the corresponding equations from the matrix. | |
| void | setValues (const labelUList &cellLabels, const UList< Type > &values) |
| Set solution in given cells to the specified values and eliminate the corresponding equations from the matrix. | |
| void | setValues (const labelUList &cellLabels, const UIndirectList< Type > &values) |
| Set solution in given cells to the specified values and eliminate the corresponding equations from the matrix. | |
| void | setReference (const label celli, const Type &value, const bool forceReference=false) |
| Set reference level for solution. | |
| void | setReferences (const labelUList &cellLabels, const Type &value, const bool forceReference=false) |
| Set reference level for solution. | |
| void | setReferences (const labelUList &cellLabels, const UList< Type > &values, const bool forceReference=false) |
| Set reference level for solution. | |
| void | setComponentReference (const label patchi, const label facei, const direction cmpt, const scalar value) |
| Set reference level for a component of the solution on a given patch face. | |
| void | addFvMatrix (fvMatrix< Type > &matrix) |
| Add fvMatrix. | |
| void | relax (const scalar alpha) |
| Relax matrix (for steady-state solution). | |
| void | relax () |
| Relax matrix (for steady-state solution). | |
| void | boundaryManipulate (typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values) |
| Manipulate based on a boundary field. | |
| autoPtr< fvSolver > | solver (const dictionary &) |
| Construct and return the solver. | |
| autoPtr< fvSolver > | solver () |
| Construct and return the solver. | |
| SolverPerformance< Type > | solveSegregatedOrCoupled (const dictionary &) |
| Solve segregated or coupled returning the solution statistics. | |
| SolverPerformance< Type > | solveSegregatedOrCoupled () |
| Solve segregated or coupled returning the solution statistics. | |
| SolverPerformance< Type > | solveSegregated (const dictionary &) |
| Solve segregated returning the solution statistics. | |
| SolverPerformance< Type > | solveCoupled (const dictionary &) |
| Solve coupled returning the solution statistics. | |
| SolverPerformance< Type > | solve (const dictionary &) |
| Solve returning the solution statistics. | |
| SolverPerformance< Type > | solve (const word &name) |
| Solve returning the solution statistics. | |
| SolverPerformance< Type > | solve () |
| Solve returning the solution statistics. | |
| tmp< Field< Type > > | residual () const |
| Return the matrix residual. | |
| tmp< scalarField > | D () const |
| Return the matrix scalar diagonal. | |
| tmp< Field< Type > > | DD () const |
| Return the matrix Type diagonal. | |
| tmp< volScalarField > | A () const |
| Return the central coefficient. | |
| tmp< GeometricField< Type, fvPatchField, volMesh > > | H () const |
| Return the H operation source. | |
| tmp< volScalarField > | H1 () const |
| Return H(1). | |
| tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > | flux () const |
| Return the face-flux field from the matrix. | |
| const dictionary & | solverDict (const word &name) const |
Return the solver dictionary (from fvSolution) for name. | |
| const dictionary & | solverDict () const |
| Return the solver dictionary for psi, taking into account finalIteration. | |
| void | operator= (const fvMatrix< Type > &) |
| void | operator= (const tmp< fvMatrix< Type > > &) |
| void | negate () |
| Inplace negate. | |
| void | operator+= (const fvMatrix< Type > &) |
| void | operator+= (const tmp< fvMatrix< Type > > &) |
| void | operator-= (const fvMatrix< Type > &) |
| void | operator-= (const tmp< fvMatrix< Type > > &) |
| void | operator+= (const DimensionedField< Type, volMesh > &) |
| void | operator+= (const tmp< DimensionedField< Type, volMesh > > &) |
| void | operator+= (const tmp< GeometricField< Type, fvPatchField, volMesh > > &) |
| void | operator-= (const DimensionedField< Type, volMesh > &) |
| void | operator-= (const tmp< DimensionedField< Type, volMesh > > &) |
| void | operator-= (const tmp< GeometricField< Type, fvPatchField, volMesh > > &) |
| void | operator+= (const dimensioned< Type > &) |
| void | operator-= (const dimensioned< Type > &) |
| void | operator+= (Foam::zero) |
| void | operator-= (Foam::zero) |
| void | operator*= (const volScalarField::Internal &) |
| void | operator*= (const tmp< volScalarField::Internal > &) |
| void | operator*= (const tmp< volScalarField > &) |
| void | operator*= (const dimensioned< scalar > &) |
| template<typename E> | |
| fvMatrix (const GeometricField< Type, fvPatchField, volMesh > &psi, const Expression::fvMatrixExpression< E, typename E::DiagExpr, typename E::UpperExpr, typename E::LowerExpr, typename E::FaceFluxExpr, typename E::SourceExpr > &expr) | |
| Construct given a field to solve for and expressions for all the components. | |
| Expression::fvMatrixConstRefWrap< fvMatrix< Type > > | expr () const |
| Wrap value as expression. | |
| template<typename E> | |
| void | operator= (const Expression::fvMatrixExpression< E, typename E::DiagExpr, typename E::UpperExpr, typename E::LowerExpr, typename E::FaceFluxExpr, typename E::SourceExpr > &expr) |
| Assign values from expression. | |
| void | setComponentReference (const label patchi, const label facei, const direction, const scalar value) |
| Foam::autoPtr< Foam::fvMatrix< Foam::scalar >::fvSolver > | solver (const dictionary &solverControls) |
| Foam::solverPerformance | solveSegregated (const dictionary &solverControls) |
| Foam::tmp< Foam::scalarField > | residual () const |
| Foam::tmp< Foam::volScalarField > | H () const |
| Foam::tmp< Foam::volScalarField > | H1 () const |
| void | setComponentReference (const label patchi, const label facei, const direction, const scalar value) |
| autoPtr< fvMatrix< scalar >::fvSolver > | solver (const dictionary &) |
| solverPerformance | solveSegregated (const dictionary &) |
| tmp< scalarField > | residual () const |
| tmp< volScalarField > | H () const |
| tmp< volScalarField > | H1 () const |
| Public Member Functions inherited from refCount | |
| constexpr | refCount () noexcept |
| Default construct, initializing count to 0. | |
| int | use_count () const noexcept |
| Return the current reference count. | |
| bool | unique () const noexcept |
| Return true if the reference count is zero. | |
| void | operator++ () noexcept |
| Increment the reference count. | |
| void | operator++ (int) noexcept |
| Increment the reference count. | |
| void | operator-- () noexcept |
| Decrement the reference count. | |
| void | operator-- (int) noexcept |
| Decrement the reference count. | |
| Public Member Functions inherited from lduMatrix | |
| ClassName ("lduMatrix") | |
| lduMatrix (const lduMesh &mesh) | |
| Construct (without coefficients) for an LDU addressed mesh. | |
| lduMatrix (const lduMatrix &) | |
| Copy construct. | |
| lduMatrix (lduMatrix &&) | |
| Move construct. | |
| lduMatrix (lduMatrix &, bool reuse) | |
| Construct as copy or re-use as specified. | |
| lduMatrix (const lduMesh &mesh, Istream &is) | |
| Construct given an LDU addressed mesh and an Istream from which the coefficients are read. | |
| ~lduMatrix ()=default | |
| Destructor. | |
| const lduMesh & | mesh () const noexcept |
| Return the LDU mesh from which the addressing is obtained. | |
| void | setLduMesh (const lduMesh &m) |
| Set the LDU mesh containing the addressing. | |
| const lduAddressing & | lduAddr () const |
| Return the LDU addressing. | |
| const lduSchedule & | patchSchedule () const |
| Return the patch evaluation schedule. | |
| 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) |
| bool | hasLowerCSR () const noexcept |
| const scalarField & | lowerCSR () const |
| scalarField & | lowerCSR () |
| const solveScalarField & | work () const |
| Work array. | |
| solveScalarField & | work (label size) const |
| Work array. | |
| word | matrixTypeName () const |
| The matrix type (empty, diagonal, symmetric, ...). | |
| bool | hasDiag () const noexcept |
| bool | hasUpper () const noexcept |
| bool | hasLower () const noexcept |
| bool | diagonal () const noexcept |
| Matrix has diagonal only. | |
| bool | symmetric () const noexcept |
| Matrix is symmetric. | |
| bool | asymmetric () const noexcept |
| Matrix is asymmetric (ie, full). | |
| 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. | |
| void | Tmul (solveScalarField &, const tmp< solveScalarField > &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &, const direction cmpt) const |
| Matrix transpose multiplication with updated interfaces. | |
| void | sumA (solveScalarField &, const FieldField< Field, scalar > &, const lduInterfaceFieldPtrsList &) const |
| Sum the coefficients on each row of the matrix. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| void | operator= (const lduMatrix &) |
| Copy assignment. | |
| void | operator= (lduMatrix &&) |
| Move assignment. | |
| 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 |
Protected Member Functions | |
| template<class Type2> | |
| void | addToInternalField (const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const |
| Add patch contribution to internal field. | |
| template<class Type2> | |
| void | addToInternalField (const labelUList &addr, const tmp< Field< Type2 > > &tpf, Field< Type2 > &intf) const |
| template<class Type2> | |
| void | subtractFromInternalField (const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const |
| Subtract patch contribution from internal field. | |
| template<class Type2> | |
| void | subtractFromInternalField (const labelUList &addr, const tmp< Field< Type2 > > &tpf, Field< Type2 > &intf) const |
| bool | checkImplicit (const label fieldi=0) |
| Name the implicit assembly addressing. | |
| void | addBoundaryDiag (scalarField &diag, const direction cmpt) const |
| void | addCmptAvBoundaryDiag (scalarField &diag) const |
| void | addBoundarySource (Field< Type > &source, const bool couples=true) const |
| template<template< class > class ListType> | |
| void | setValuesFromList (const labelUList &cellLabels, const ListType< Type > &values) |
| Set solution in given cells to the specified values. | |
Friends | |
| class | fvSolver |
| Declare friendship with the fvSolver class. | |
| tmp< GeometricField< Type, fvPatchField, volMesh > > | operator& (const fvMatrix< Type > &, const DimensionedField< Type, volMesh > &) |
| tmp< GeometricField< Type, fvPatchField, volMesh > > | operator& (const fvMatrix< Type > &, const tmp< GeometricField< Type, fvPatchField, volMesh > > &) |
| tmp< GeometricField< Type, fvPatchField, volMesh > > | operator& (const tmp< fvMatrix< Type > > &, const DimensionedField< Type, volMesh > &) |
| tmp< GeometricField< Type, fvPatchField, volMesh > > | operator& (const tmp< fvMatrix< Type > > &, const tmp< GeometricField< Type, fvPatchField, volMesh > > &) |
| Ostream & | operator<< (Ostream &, const fvMatrix< Type > &) |
Additional Inherited Members | |
| Static Public Attributes inherited from lduMatrix | |
| static const Enum< normTypes > | normTypesNames_ |
| Names for the normTypes. | |
| static constexpr const label | defaultMaxIter = 1000 |
| Default maximum number of iterations for solvers (1000). | |
| static const scalar | defaultTolerance = 1e-6 |
| Default (absolute) tolerance (1e-6). | |
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition at line 114 of file fvMatrix.H.
| typedef GeometricField<Type, fvPatchField, volMesh> psiFieldType |
The geometric field type for psi.
Definition at line 128 of file fvMatrix.H.
| typedef GeometricField<Type, fvsPatchField, surfaceMesh> faceFluxFieldType |
Field type for face flux (for non-orthogonal correction).
Definition at line 135 of file fvMatrix.H.
| typedef std::unique_ptr<faceFluxFieldType> faceFluxFieldPtrType |
Declare return type of the faceFluxCorrectionPtr() function.
Definition at line 584 of file fvMatrix.H.
| fvMatrix | ( | const GeometricField< Type, fvPatchField, volMesh > & | psi, |
| const dimensionSet & | ds ) |
Construct given a field to solve for.
Definition at line 350 of file fvMatrix.C.
References boundary, checkImplicit(), DebugInFunction, Foam::endl(), forAll, lduMatrix::lduMatrix(), lduMatrix::mesh(), psi(), psi, and Foam::Zero.
Referenced by addFvMatrix(), fvMatrix(), fvMatrix(), fvMatrix< Type >::fvSolver::fvSolver(), operator+=(), operator+=(), operator-=(), operator-=(), operator=(), and operator=().


| fvMatrix | ( | const fvMatrix< Type > & | fvm | ) |
Copy construct.
Definition at line 394 of file fvMatrix.C.
References DebugInFunction, Foam::endl(), fvMatrix(), and lduMatrix::lduMatrix().

| fvMatrix | ( | const tmp< fvMatrix< Type > > & | tmat | ) |
Copy/move construct from tmp<fvMatrix<Type>>.
Definition at line 420 of file fvMatrix.C.
References DebugInFunction, Foam::endl(), fvMatrix(), and lduMatrix::lduMatrix().

|
inline |
Deprecated(2022-05) - construct with dimensionSet instead.
Definition at line 364 of file fvMatrix.H.
|
virtual |
Destructor.
Definition at line 458 of file fvMatrix.C.
References DebugInFunction, and Foam::endl().

| fvMatrix | ( | const GeometricField< Type, fvPatchField, volMesh > & | psi, |
| const Expression::fvMatrixExpression< E, typename E::DiagExpr, typename E::UpperExpr, typename E::LowerExpr, typename E::FaceFluxExpr, typename E::SourceExpr > & | expr ) |
Construct given a field to solve for and expressions for all the components.
Definition at line 2984 of file fvMatrix.C.
References boundary, checkImplicit(), DebugInFunction, dimensions(), Foam::endl(), expr(), lduMatrix::lduMatrix(), lduMatrix::mesh(), psi(), psi, and Foam::Zero.

|
protected |
Add patch contribution to internal field.
Definition at line 40 of file fvMatrix.C.
References Foam::abort(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, forAll, and UList< T >::size().
Referenced by addBoundaryDiag(), addBoundarySource(), addCmptAvBoundaryDiag(), addToInternalField(), DD(), and H1().


|
protected |
Definition at line 64 of file fvMatrix.C.
References addToInternalField().

|
protected |
Subtract patch contribution from internal field.
Definition at line 78 of file fvMatrix.C.
References Foam::abort(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, forAll, and UList< T >::size().
Referenced by subtractFromInternalField().


|
protected |
Definition at line 102 of file fvMatrix.C.
References subtractFromInternalField().

|
protected |
Name the implicit assembly addressing.
Definition at line 314 of file fvMatrix.C.
References Foam::endl(), forAll, Foam::name(), Foam::Pout, and psi.
Referenced by addFvMatrix(), fvMatrix(), and fvMatrix().


|
protected |
Definition at line 115 of file fvMatrix.C.
References addToInternalField(), Foam::component(), lduMatrix::diag(), forAll, globalPatchID(), lduMatrix::lduAddr(), nMatrices(), and psi.
Referenced by H(), residual(), fvMatrix< scalar >::residual(), solveCoupled(), solveSegregated(), and fvMatrix< scalar >::solveSegregated().


|
protected |
Definition at line 144 of file fvMatrix.C.
References addToInternalField(), Foam::cmptAv(), lduMatrix::diag(), forAll, globalPatchID(), lduMatrix::lduAddr(), nMatrices(), and psi.


|
protected |
Definition at line 168 of file fvMatrix.C.
References addToInternalField(), Foam::cmptMultiply(), fvPatchFieldBase::coupled(), forAll, globalPatchID(), lduMatrix::lduAddr(), nMatrices(), fvPatchField< Type >::patchNeighbourField(), psi, and source().
Referenced by H(), fvMatrix< scalar >::H(), residual(), fvMatrix< scalar >::residual(), solveCoupled(), solveSegregated(), and fvMatrix< scalar >::solveSegregated().


|
protected |
Set solution in given cells to the specified values.
Definition at line 218 of file fvMatrix.C.
References lduMatrix::asymmetric(), cells, lduMatrix::diag(), forAll, lduMatrix::lower(), lduMatrix::mesh(), primitiveFieldRef(), psi(), lduMatrix::symmetric(), lduMatrix::upper(), and Foam::Zero.
Referenced by setValues(), setValues(), and setValues().


| ClassName | ( | "fvMatrix< Type >" | ) |
Construct and return a clone.
Definition at line 376 of file fvMatrix.H.
Referenced by setBounAndInterCoeffs().

|
inline |
Definition at line 394 of file fvMatrix.H.
Referenced by addBoundaryDiag(), addBoundarySource(), addCmptAvBoundaryDiag(), addFvMatrix(), createOrUpdateLduPrimitiveAssembly(), flux(), manipulateMatrix(), setBounAndInterCoeffs(), setInterfaces(), fvMatrix< scalar >::solveSegregated(), and transferFvMatrixCoeffs().

|
inline |
Definition at line 399 of file fvMatrix.H.
Referenced by addFvMatrix(), fvMatrix< scalar >::psi(), fvMatrix< scalar >::psi(), setBounAndInterCoeffs(), and transferFvMatrixCoeffs().

|
inline |
Definition at line 404 of file fvMatrix.H.
|
inline |
Definition at line 409 of file fvMatrix.H.
Referenced by addBoundaryDiag(), addBoundarySource(), addCmptAvBoundaryDiag(), and setInterfaces().

| void transferFvMatrixCoeffs | ( | ) |
Transfer lower, upper, diag and source to this fvMatrix.
Definition at line 814 of file fvMatrix.C.
References lduMatrix::asymmetric(), Foam::diag(), Foam::faceMap(), forAll, lduMeshPtr(), lduMatrix::lower(), matrix(), nMatrices(), List< T >::setSize(), source(), lduMatrix::upper(), and Foam::Zero.
Referenced by fvMatrix< scalar >::solveSegregated().


| void createOrUpdateLduPrimitiveAssembly | ( | ) |
Create or update ldu assembly.
Definition at line 907 of file fvMatrix.C.
References lduAddressing::clearOut(), Foam::endl(), Foam::Info, io, lduPrimitiveMesh::lduAddr(), lduMeshPtr(), DimensionedField< Type, GeoMesh >::mesh(), nMatrices(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, psi, IOobjectOption::REGISTER, UPtrList< T >::set(), regIOobject::store(), and lduPrimitiveMeshAssembly::update().
Referenced by fvMatrix< scalar >::solveSegregated().


| Foam::lduPrimitiveMeshAssembly * lduMeshPtr | ( | ) |
Access to lduPrimitiveMeshAssembly.
Definition at line 881 of file fvMatrix.C.
Referenced by createOrUpdateLduPrimitiveAssembly(), fvMatrix< scalar >::globalPatchID(), fvMatrix< scalar >::lduMeshAssembly(), manipulateMatrix(), mapContributions(), setBounAndInterCoeffs(), setInterfaces(), fvMatrix< scalar >::solveSegregated(), and transferFvMatrixCoeffs().

| const Foam::lduPrimitiveMeshAssembly * lduMeshPtr | ( | ) | const |
Const Access to lduPrimitiveMeshAssembly.
Definition at line 894 of file fvMatrix.C.
| void manipulateMatrix | ( | direction | cmp | ) |
Manipulate matrix.
Definition at line 791 of file fvMatrix.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), forAll, lduMeshPtr(), nMatrices(), and psi.
Referenced by fvMatrix< scalar >::solveSegregated().


| void setBounAndInterCoeffs | ( | ) |
Manipulate boundary/internal coeffs for coupling.
Definition at line 661 of file fvMatrix.C.
References boundary, boundaryCoeffs(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), clone(), forAll, internalCoeffs(), lduMeshPtr(), matrix(), nMatrices(), psi(), psi, List< T >::set(), List< T >::setSize(), and Foam::Zero.
Referenced by fvMatrix< scalar >::solveSegregated().


| void setInterfaces | ( | lduInterfaceFieldPtrsList & | interfaces, |
| PtrDynList< lduInterfaceField > & | newInterfaces ) |
Set interfaces.
Definition at line 470 of file fvMatrix.C.
References PtrDynList< T, SizeMin >::append(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), forAll, globalPatchID(), Foam::isA(), UPtrList< T >::last(), lduMeshPtr(), nMatrices(), psi, Foam::refCast(), UPtrList< T >::set(), and UPtrList< T >::setSize().
Referenced by fvMatrix< scalar >::solveSegregated().


| void mapContributions | ( | label | fieldi, |
| const FieldField< Field, Type > & | fluxContrib, | ||
| FieldField< Field, Type > & | contrib, | ||
| bool | internal ) const |
Add internal and boundary contribution to local patches.
Definition at line 548 of file fvMatrix.C.
References lduPrimitiveMeshAssembly::cellBoundMap(), cellId, Foam::cmptMultiply(), coupled, faceId(), lduPrimitiveMeshAssembly::facePatchFaceMap(), forAll, lduMeshPtr(), DimensionedField< Type, GeoMesh >::mesh(), lduPrimitiveMeshAssembly::patchLocalToGlobalMap(), lduPrimitiveMeshAssembly::patchMap(), pp(), psi(), psi, List< T >::setSize(), and PtrList< T >::setSize().
Referenced by flux().


|
inline |
Return optional lduAdressing.
Definition at line 478 of file fvMatrix.H.
Referenced by cyclicACMIFvPatchField< Type >::manipulateMatrix(), cyclicAMIFvPatchField< Type >::manipulateMatrix(), cyclicFvPatchField< Type >::manipulateMatrix(), and mixedEnergyFvPatchScalarField::manipulateMatrix().

|
inline |
Return psi.
Definition at line 486 of file fvMatrix.H.
Referenced by velocityDampingConstraint::addDamping(), oversetFvMeshBase::addInterpolation(), isotropic::addRegularisationTerm(), objectiveIncompressible::addSource(), objectiveNutSqr::addSource(), objectivePowerDissipation::addSource(), acousticDampingSource::addSup(), acousticDampingSource::addSup(), acousticDampingSource::addSup(), atmCoriolisUSource::addSup(), atmCoriolisUSource::addSup(), atmCoriolisUSource::addSup(), atmPlantCanopyUSource::addSup(), atmPlantCanopyUSource::addSup(), atmPlantCanopyUSource::addSup(), explicitPorositySource::addSup(), explicitPorositySource::addSup(), explicitPorositySource::addSup(), interRegionExplicitPorositySource::addSup(), interRegionExplicitPorositySource::addSup(), interRegionHeatTransferModel::addSup(), multiphaseMangrovesSource::addSup(), multiphaseMangrovesSource::addSup(), multiphaseMangrovesTurbulenceModel::addSup(), multiphaseMangrovesTurbulenceModel::addSup(), multiphaseStabilizedTurbulence::addSup(), patchCellsSource::addSup(), PhaseLimitStabilization< Type >::addSup(), radialActuationDiskSource::addSup(), radialActuationDiskSource::addSup(), radiation::addSup(), rotorDiskSource::addSup(), rotorDiskSource::addSup(), SemiImplicitSource< Type >::addSup(), solidificationMeltingSource::addSup(), topOSource::addSup(), topOSource::addSup(), topOSource::addSup(), topOSource::addSup(), Foam::checkMethod(), Foam::checkMethod(), Foam::checkMethod(), optionList::constrain(), oversetFvPatchField< Type >::fringeFlux(), fvMatrix(), fvMatrix(), cyclicACMIFvPatchField< Type >::manipulateMatrix(), cyclicAMIFvPatchField< Type >::manipulateMatrix(), mapContributions(), MomentumTransferPhaseSystem< BasePhaseSystem >::momentumTransfer(), oversetFvMeshBase::normalisation(), Foam::oversetPatchPhiErr(), setBounAndInterCoeffs(), regularisationPDE::setValues(), setValuesFromList(), fvMatrix< Type >::fvSolver::solve(), solveCoupled(), oversetFvMeshBase::solveOverset(), solveSegregated(), Foam::Expression::Sp(), Foam::Expression::Su(), Foam::Expression::SuSp(), and oversetFvMeshBase::write().
|
inline |
Definition at line 498 of file fvMatrix.H.
|
inline |
Clear multiple fvMatrices.
Definition at line 523 of file fvMatrix.H.
|
inlinenoexcept |
Definition at line 530 of file fvMatrix.H.
Referenced by directionalPressureGradientExplicitSource::addSup(), explicitPorositySource::addSup(), explicitPorositySource::addSup(), explicitPorositySource::addSup(), fanMomentumSource::addSup(), fanMomentumSource::addSup(), interRegionExplicitPorositySource::addSup(), interRegionExplicitPorositySource::addSup(), meanVelocityForce::addSup(), patchCellsSource::addSup(), rotorDiskSource::addSup(), rotorDiskSource::addSup(), SemiImplicitSource< Type >::addSup(), Foam::checkMethod(), Foam::checkMethod(), Foam::checkMethod(), flux(), and fvMatrix().

|
inlinenoexcept |
Definition at line 535 of file fvMatrix.H.
Referenced by addBoundarySource(), oversetFvMeshBase::addInterpolation(), heatExchangerSource::addSup(), interRegionExplicitPorositySource::addSup(), interRegionExplicitPorositySource::addSup(), radialActuationDiskSource::addSup(), radialActuationDiskSource::addSup(), solidificationMeltingSource::addSup(), adjointkOmegaSST::addWallFunctionTerms(), adjointOutletVelocityFluxFvPatchVectorField::manipulateMatrix(), adjointWallVelocityFvPatchVectorField::manipulateMatrix(), kaqRWallFunctionFvPatchScalarField::manipulateMatrix(), mixedEnergyFvPatchScalarField::manipulateMatrix(), operator+=(), operator+=(), operator-=(), operator-=(), relax(), setReference(), setReferences(), setReferences(), solveCoupled(), oversetFvMeshBase::solveOverset(), solveSegregated(), Foam::Expression::Su(), Foam::Expression::SuSp(), transferFvMatrixCoeffs(), and oversetFvMeshBase::write().

|
inlinenoexcept |
Definition at line 540 of file fvMatrix.H.
|
inlinenoexcept |
fvBoundary scalar field containing pseudo-matrix coeffs for internal cells
Definition at line 549 of file fvMatrix.H.
Referenced by oversetFvMeshBase::addInterpolation(), cyclicACMIFvPatchField< Type >::manipulateMatrix(), cyclicAMIFvPatchField< Type >::manipulateMatrix(), cyclicFvPatchField< Type >::manipulateMatrix(), mixedEnergyFvPatchScalarField::manipulateMatrix(), oversetFvPatchField< Type >::manipulateMatrix(), waWallFunctionFvPatchScalarField::manipulateMatrix(), oversetFvMeshBase::normalisation(), setBounAndInterCoeffs(), solveCoupled(), oversetFvMeshBase::solveOverset(), oversetFvPatchField< Type >::storeFringeCoefficients(), and oversetFvMeshBase::write().

|
inlinenoexcept |
fvBoundary scalar field containing pseudo-matrix coeffs for internal cells
Definition at line 558 of file fvMatrix.H.
|
inlinenoexcept |
fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells
Definition at line 567 of file fvMatrix.H.
Referenced by oversetFvMeshBase::addInterpolation(), cyclicACMIFvPatchField< Type >::manipulateMatrix(), cyclicAMIFvPatchField< Type >::manipulateMatrix(), cyclicFvPatchField< Type >::manipulateMatrix(), oversetFvPatchField< Type >::manipulateMatrix(), waWallFunctionFvPatchScalarField::manipulateMatrix(), setBounAndInterCoeffs(), solveCoupled(), oversetFvMeshBase::solveOverset(), oversetFvPatchField< Type >::storeFringeCoefficients(), and oversetFvMeshBase::write().

|
inlinenoexcept |
fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells
Definition at line 576 of file fvMatrix.H.
|
inline |
Return pointer to face-flux non-orthogonal correction field.
Definition at line 589 of file fvMatrix.H.
|
inline |
Set pointer to face-flux non-orthogonal correction field.
Definition at line 597 of file fvMatrix.H.
|
inlinenoexcept |
True if face-flux non-orthogonal correction field exists.
Definition at line 605 of file fvMatrix.H.
| void setValues | ( | const labelUList & | cellLabels, |
| const Type & | value ) |
Set solution in given cells to the specified value and eliminate the corresponding equations from the matrix.
Definition at line 970 of file fvMatrix.C.
References setValuesFromList().
Referenced by fixedTemperatureConstraint::constrain(), FixedValueConstraint< Type >::constrain(), MapFieldConstraint< Type >::constrain(), epsilonWallFunctionFvPatchScalarField::manipulateMatrix(), fixedInternalValueFvPatchField< Type >::manipulateMatrix(), omegaWallFunctionFvPatchScalarField::manipulateMatrix(), and regularisationPDE::setValues().


| void setValues | ( | const labelUList & | cellLabels, |
| const UList< Type > & | values ) |
Set solution in given cells to the specified values and eliminate the corresponding equations from the matrix.
Definition at line 981 of file fvMatrix.C.
References setValuesFromList().

| void setValues | ( | const labelUList & | cellLabels, |
| const UIndirectList< Type > & | values ) |
Set solution in given cells to the specified values and eliminate the corresponding equations from the matrix.
Definition at line 992 of file fvMatrix.C.
References setValuesFromList().

| void setReference | ( | const label | celli, |
| const Type & | value, | ||
| const bool | forceReference = false ) |
Set reference level for solution.
Definition at line 1003 of file fvMatrix.C.
References Foam::diag(), and source().
Referenced by Foam::CorrectPhi(), simple::mainIter(), and while().


| void setReferences | ( | const labelUList & | cellLabels, |
| const Type & | value, | ||
| const bool | forceReference = false ) |
Set reference level for solution.
Definition at line 1019 of file fvMatrix.C.
References cellId, Foam::diag(), forAll, and source().

| void setReferences | ( | const labelUList & | cellLabels, |
| const UList< Type > & | values, | ||
| const bool | forceReference = false ) |
Set reference level for solution.
Definition at line 1042 of file fvMatrix.C.
References cellId, Foam::diag(), forAll, and source().

| void setComponentReference | ( | const label | patchi, |
| const label | facei, | ||
| const direction | cmpt, | ||
| const scalar | value ) |
Set reference level for a component of the solution on a given patch face.
Definition at line 30 of file fvMatrixSolve.C.
References lduMatrix::diag(), and UPstream::master().

| void addFvMatrix | ( | fvMatrix< Type > & | matrix | ) |
Add fvMatrix.
Definition at line 1065 of file fvMatrix.C.
References Foam::abort(), checkImplicit(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, fvMatrix(), matrix(), and nMatrices().

| void relax | ( | const scalar | alpha | ) |
Relax matrix (for steady-state solution).
alpha = 1 : diagonally equal alpha < 1 : diagonally dominant alpha = 0 : do nothing Note: Requires positive diagonal.
Definition at line 1094 of file fvMatrix.C.
References alpha, Foam::cmptMag(), Foam::cmptMax(), Foam::cmptMin(), Foam::component(), fvPatchFieldBase::coupled(), D(), DebugInFunction, Foam::diag(), Foam::endl(), forAll, InfoInFunction, lduMatrix::lduAddr(), Foam::mag(), Foam::max(), UPstream::msgType(), Foam::nl, Foam::reduce(), Foam::returnReduce(), UList< T >::size(), source(), lduMatrix::sumMagOffDiag(), and Foam::Zero.
Referenced by hydrostaticPressure::calculateAndWrite(), IATE::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), advectionDiffusion::correct(), radiativeIntensityRay::correct(), kineticTheoryModel::correct(), thixotropicViscosity::correct(), waxSolventEvaporation::correctModel(), thermo::evolveRegion(), age::execute(), electricPotential::execute(), energyTransport::execute(), scalarTransport::execute(), adjointSimple::mainIter(), adjointEikonalSolver::solve(), populationBalanceModel::solve(), twoPhaseSystem::solve(), reactingOneDim::solveEnergy(), and thermalBaffle::solveEnergy().


| void relax | ( | ) |
Relax matrix (for steady-state solution).
alpha is read from controlDict
Definition at line 1242 of file fvMatrix.C.
References Foam::name(), and relax().

| void boundaryManipulate | ( | typename GeometricField< Type, fvPatchField, volMesh >::Boundary & | values | ) |
Manipulate based on a boundary field.
Definition at line 1259 of file fvMatrix.C.
References forAll.
Referenced by adjointSimple::mainIter(), adjointMeshMovementSolver::solve(), multiphaseSystem::solveAlphas(), and oversetFvMeshBase::solveOverset().

| autoPtr< fvSolver > solver | ( | const dictionary & | ) |
Construct and return the solver.
Use the given solver controls
Referenced by radiativeIntensityRay::correct(), and solver().

| Foam::autoPtr< typename Foam::fvMatrix< Type >::fvSolver > solver | ( | ) |
Construct and return the solver.
Solver controls read from fvSolution
Definition at line 327 of file fvMatrixSolve.C.
References solver(), and solverDict().

| Foam::SolverPerformance< Type > solveSegregatedOrCoupled | ( | const dictionary & | solverControls | ) |
Solve segregated or coupled returning the solution statistics.
Use the given solver controls
Definition at line 54 of file fvMatrixSolve.C.
References addProfiling, polyMesh::defaultRegion, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, dictionary::getOrDefault(), Foam::Info, mesh, regionName, solve(), solveCoupled(), and solveSegregated().
Referenced by fvMesh::solve(), fvMesh::solve(), fvMesh::solve(), fvMesh::solve(), fvMesh::solve(), surfaceAlignedSBRStressFvMotionSolver::solve(), and solveSegregatedOrCoupled().


| Foam::SolverPerformance< Type > solveSegregatedOrCoupled | ( | ) |
Solve segregated or coupled returning the solution statistics.
Solver controls read from fvSolution
Definition at line 309 of file fvMatrixSolve.C.
References solverDict(), and solveSegregatedOrCoupled().

| Foam::SolverPerformance< Type > solveSegregated | ( | const dictionary & | solverControls | ) |
Solve segregated returning the solution statistics.
Use the given solver controls
Definition at line 104 of file fvMatrixSolve.C.
References addBoundaryDiag(), addBoundarySource(), lduMatrix::diag(), Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, dictionary::getOrDefault(), Foam::Info, lduMatrix::initMatrixInterfaces(), mesh, lduMatrix::solver::New(), UPstream::nRequests(), SolverPerformance< Type >::print(), psi(), refPtr< T >::ref(), SolverPerformance< Type >::replace(), SolverPerformance< Type >::solverName(), source(), and lduMatrix::updateMatrixInterfaces().
Referenced by solveSegregatedOrCoupled().


| Foam::SolverPerformance< Type > solveCoupled | ( | const dictionary & | solverControls | ) |
Solve coupled returning the solution statistics.
Use the given solver controls
Definition at line 243 of file fvMatrixSolve.C.
References addBoundaryDiag(), addBoundarySource(), boundaryCoeffs(), LduMatrix< Type, DType, LUType >::diag(), lduMatrix::diag(), Foam::endl(), dictionary::getOrDefault(), Foam::Info, LduMatrix< Type, DType, LUType >::interfaces(), LduMatrix< Type, DType, LUType >::interfacesLower(), LduMatrix< Type, DType, LUType >::interfacesUpper(), internalCoeffs(), LduMatrix< Type, DType, LUType >::lower(), lduMatrix::lower(), mesh, LduMatrix< Type, DType, LUType >::solver::New(), SolverPerformance< Type >::print(), psi(), source(), LduMatrix< Type, DType, LUType >::source(), LduMatrix< Type, DType, LUType >::upper(), and lduMatrix::upper().
Referenced by solveSegregatedOrCoupled().


| Foam::SolverPerformance< Type > solve | ( | const dictionary & | solverControls | ) |
Solve returning the solution statistics.
Use the given solver controls
Definition at line 316 of file fvMatrixSolve.C.
Referenced by Implicit< CloudType >::cacheFields(), hydrostaticPressure::calculateAndWrite(), IATE::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), advectionDiffusion::correct(), kineticTheoryModel::correct(), thixotropicViscosity::correct(), waxSolventEvaporation::correctModel(), Foam::CorrectPhi(), Foam::CorrectPhi(), thermo::evolveRegion(), age::execute(), electricPotential::execute(), energyTransport::execute(), scalarTransport::execute(), simple::mainIter(), adjointEikonalSolver::solve(), populationBalanceModel::solve(), elasticityMotionSolver::solve(), laplacianMotionSolver::solve(), pLaplacianMotionSolver::solve(), twoPhaseSystem::solve(), multiphaseSystem::solveAlphas(), reactingOneDim::solveContinuity(), reactingOneDim::solveEnergy(), thermalBaffle::solveEnergy(), Helmholtz::solveEqn(), shapeDesignVariables::solveMeshMovementEqn(), solveSegregatedOrCoupled(), reactingOneDim::solveSpeciesMass(), kinematicSingleLayer::solveThickness(), MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::solveYi(), and Foam::fvc::spreadSource().

| Foam::SolverPerformance< Type > solve | ( | const word & | name | ) |
Solve returning the solution statistics.
Uses name solver controls from fvSolution
Definition at line 341 of file fvMatrixSolve.C.
References Foam::name(), solve(), and solverDict().

| Foam::SolverPerformance< Type > solve | ( | ) |
Solve returning the solution statistics.
Solver controls read from fvSolution
Definition at line 348 of file fvMatrixSolve.C.
References solve(), and solverDict().
Referenced by solve(), and solve().


| Foam::tmp< Foam::Field< Type > > residual | ( | ) | const |
Return the matrix residual.
Definition at line 355 of file fvMatrixSolve.C.
References addBoundaryDiag(), addBoundarySource(), Field< Type >::component(), Foam::New(), lduMatrix::residual(), and Foam::Zero.

| Foam::tmp< Foam::scalarField > D | ( | ) | const |
Return the matrix scalar diagonal.
Definition at line 1273 of file fvMatrix.C.
References addCmptAvBoundaryDiag(), Foam::diag(), and tmp< T >::New().
Referenced by relax().


| Foam::tmp< Foam::Field< Type > > DD | ( | ) | const |
Return the matrix Type diagonal.
Definition at line 1282 of file fvMatrix.C.
References addToInternalField(), fvPatchFieldBase::coupled(), Foam::diag(), forAll, lduMatrix::lduAddr(), tmp< T >::ref(), and UList< T >::size().

| Foam::tmp< Foam::volScalarField > A | ( | ) | const |
Return the central coefficient.
Definition at line 1306 of file fvMatrix.C.
References D, Foam::dimVol, fvPatchFieldBase::extrapolatedCalculatedType(), and GeometricField< scalar, fvPatchField, volMesh >::New().
Referenced by directionalPressureGradientExplicitSource::constrain(), meanVelocityForce::constrain(), and adjointSimple::mainIter().


| Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > H | ( | ) | const |
Return the H operation source.
Definition at line 1325 of file fvMatrix.C.
References addBoundaryDiag(), addBoundarySource(), addCmptAvBoundaryDiag(), Foam::dimVol, fvPatchFieldBase::extrapolatedCalculatedType(), lduMatrix::H(), Field< Type >::negate(), GeometricField< Type, PatchField, GeoMesh >::New(), and Foam::Zero.
Referenced by adjointSimple::mainIter().


| Foam::tmp< Foam::volScalarField > H1 | ( | ) | const |
Return H(1).
Definition at line 1377 of file fvMatrix.C.
References addToInternalField(), Foam::component(), fvPatchFieldBase::coupled(), Foam::dimVol, fvPatchFieldBase::extrapolatedCalculatedType(), forAll, lduMatrix::H1(), lduMatrix::lduAddr(), GeometricField< scalar, fvPatchField, volMesh >::New(), and UList< T >::size().
Referenced by adjointSimple::mainIter().


| Foam::tmp< Foam::GeometricField< Type, Foam::fvsPatchField, Foam::surfaceMesh > > flux | ( | ) | const |
Return the face-flux field from the matrix.
Definition at line 1414 of file fvMatrix.C.
References Foam::abort(), Foam::cmptMultiply(), dimensions(), lduMatrix::faceH(), Foam::FatalError, FatalErrorInFunction, forAll, mapContributions(), GeometricField< Type, PatchField, GeoMesh >::New(), and nMatrices().
Referenced by Implicit< CloudType >::cacheFields(), incompressiblePrimalSolver::correctBoundaryConditions(), Foam::CorrectPhi(), Foam::CorrectPhi(), scalarTransport::execute(), simple::mainIter(), twoPhaseSystem::solve(), multiphaseSystem::solveAlphas(), and kinematicSingleLayer::solveThickness().


| const Foam::dictionary & solverDict | ( | const word & | name | ) | const |
Return the solver dictionary (from fvSolution) for name.
Definition at line 1512 of file fvMatrix.C.
References Foam::name().
Referenced by solver(), and solveSegregatedOrCoupled().


| const Foam::dictionary & solverDict | ( | ) | const |
Return the solver dictionary for psi, taking into account finalIteration.
Definition at line 1522 of file fvMatrix.C.
Referenced by solve(), and solve().

| void operator= | ( | const fvMatrix< Type > & | fvmv | ) |
Definition at line 1534 of file fvMatrix.C.
References Foam::abort(), Foam::FatalError, FatalErrorInFunction, fvMatrix(), and lduMatrix::operator=().
Referenced by operator=().


Definition at line 1572 of file fvMatrix.C.
References fvMatrix(), and operator=().

| void negate | ( | ) |
Inplace negate.
Definition at line 1580 of file fvMatrix.C.
References lduMatrix::negate().

| void operator+= | ( | const fvMatrix< Type > & | fvmv | ) |
Definition at line 1595 of file fvMatrix.C.
References Foam::checkMethod(), fvMatrix(), and lduMatrix::operator+=().
Referenced by operator+=(), operator+=(), and operator+=().


Definition at line 1624 of file fvMatrix.C.
References fvMatrix(), and operator+=().

| void operator-= | ( | const fvMatrix< Type > & | fvmv | ) |
Definition at line 1632 of file fvMatrix.C.
References Foam::checkMethod(), fvMatrix(), and lduMatrix::operator-=().
Referenced by operator-=(), operator-=(), and operator-=().


Definition at line 1661 of file fvMatrix.C.
References fvMatrix(), and operator-=().

| void operator+= | ( | const DimensionedField< Type, volMesh > & | su | ) |
Definition at line 1669 of file fvMatrix.C.
References Foam::checkMethod(), and source().

| void operator+= | ( | const tmp< DimensionedField< Type, volMesh > > & | tsu | ) |
Definition at line 1680 of file fvMatrix.C.
References operator+=().

| void operator+= | ( | const tmp< GeometricField< Type, fvPatchField, volMesh > > & | tsu | ) |
Definition at line 1691 of file fvMatrix.C.
References operator+=().

| void operator-= | ( | const DimensionedField< Type, volMesh > & | su | ) |
Definition at line 1702 of file fvMatrix.C.
References Foam::checkMethod(), and source().

| void operator-= | ( | const tmp< DimensionedField< Type, volMesh > > & | tsu | ) |
Definition at line 1713 of file fvMatrix.C.
References operator-=().

| void operator-= | ( | const tmp< GeometricField< Type, fvPatchField, volMesh > > & | tsu | ) |
Definition at line 1724 of file fvMatrix.C.
References operator-=().

| void operator+= | ( | const dimensioned< Type > & | su | ) |
Definition at line 1735 of file fvMatrix.C.
References DimensionedField< Type, GeoMesh >::mesh(), psi, and source().

| void operator-= | ( | const dimensioned< Type > & | su | ) |
Definition at line 1745 of file fvMatrix.C.
References DimensionedField< Type, GeoMesh >::mesh(), psi, and source().

|
inline |
Definition at line 860 of file fvMatrix.H.
|
inline |
Definition at line 861 of file fvMatrix.H.
| void operator*= | ( | const volScalarField::Internal & | dsf | ) |
Definition at line 1755 of file fvMatrix.C.
References Foam::abort(), Foam::FatalError, FatalErrorInFunction, forAll, and lduMatrix::operator*=().
Referenced by operator*=(), and operator*=().


| void operator*= | ( | const tmp< volScalarField::Internal > & | tfld | ) |
Definition at line 1785 of file fvMatrix.C.
References operator*=().

| void operator*= | ( | const tmp< volScalarField > & | tfld | ) |
Definition at line 1796 of file fvMatrix.C.
References operator*=().

| void operator*= | ( | const dimensioned< scalar > & | ds | ) |
Definition at line 1807 of file fvMatrix.C.
References lduMatrix::operator*=().

| Foam::Expression::fvMatrixConstRefWrap< Foam::fvMatrix< Type > > expr | ( | ) | const |
Wrap value as expression.
Definition at line 3026 of file fvMatrix.C.
Referenced by fvMatrix(), and operator=().

| void operator= | ( | const Expression::fvMatrixExpression< E, typename E::DiagExpr, typename E::UpperExpr, typename E::LowerExpr, typename E::FaceFluxExpr, typename E::SourceExpr > & | expr | ) |
Assign values from expression.
Definition at line 3034 of file fvMatrix.C.
References expr().

| void setComponentReference | ( | const label | patchi, |
| const label | facei, | ||
| const direction | , | ||
| const scalar | value ) |
Definition at line 33 of file fvScalarMatrix.C.
| Foam::autoPtr< Foam::fvMatrix< Foam::scalar >::fvSolver > solver | ( | const dictionary & | solverControls | ) |
Definition at line 58 of file fvScalarMatrix.C.
| Foam::solverPerformance solveSegregated | ( | const dictionary & | solverControls | ) |
Definition at line 154 of file fvScalarMatrix.C.
| Foam::tmp< Foam::scalarField > residual | ( | ) | const |
Definition at line 308 of file fvScalarMatrix.C.
| Foam::tmp< Foam::volScalarField > H | ( | ) | const |
Definition at line 337 of file fvScalarMatrix.C.
| Foam::tmp< Foam::volScalarField > H1 | ( | ) | const |
Definition at line 359 of file fvScalarMatrix.C.
| void setComponentReference | ( | const label | patchi, |
| const label | facei, | ||
| const direction | , | ||
| const scalar | value ) |
| autoPtr< fvMatrix< scalar >::fvSolver > solver | ( | const dictionary & | ) |
| solverPerformance solveSegregated | ( | const dictionary & | ) |
| tmp< scalarField > residual | ( | ) | const |
| tmp< volScalarField > H | ( | ) | const |
| tmp< volScalarField > H1 | ( | ) | const |
|
friend |
Declare friendship with the fvSolver class.
Definition at line 205 of file fvMatrix.H.
|
friend |
|
friend |
|
friend |
|
friend |