fvMatrix< Type > Class Template Reference

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...

Classes

class  fvSolver
 

Public Types

typedef GeometricField< Type, fvsPatchField, surfaceMesh > * surfaceTypeFieldPtr
 Declare return type of the faceFluxCorrectionPtr() function. More...
 

Public Member Functions

 ClassName ("fvMatrix")
 
 fvMatrix (const GeometricField< Type, fvPatchField, volMesh > &psi, const dimensionSet &ds)
 Construct given a field to solve for. More...
 
 fvMatrix (const fvMatrix< Type > &)
 Copy construct. More...
 
 fvMatrix (const tmp< fvMatrix< Type >> &)
 Copy/move construct from tmp<fvMatrix<Type>> More...
 
 fvMatrix (const GeometricField< Type, fvPatchField, volMesh > &psi, Istream &is)
 Construct from Istream given field to solve for. More...
 
tmp< fvMatrix< Type > > clone () const
 Clone. More...
 
virtual ~fvMatrix ()
 Destructor. More...
 
const GeometricField< Type, fvPatchField, volMesh > & psi () const
 
const dimensionSetdimensions () const
 
Field< Type > & source ()
 
const Field< Type > & source () const
 
const FieldField< Field, Type > & internalCoeffs () const
 
FieldField< Field, Type > & internalCoeffs ()
 
const FieldField< Field, Type > & boundaryCoeffs () const
 
FieldField< Field, Type > & boundaryCoeffs ()
 
surfaceTypeFieldPtrfaceFluxCorrectionPtr ()
 Return pointer to face-flux non-orthogonal correction field. More...
 
void setValues (const labelUList &cellLabels, const Type &value)
 
void setValues (const labelUList &cellLabels, const UList< Type > &values)
 
void setValues (const labelUList &cellLabels, const UIndirectList< Type > &values)
 
void setReference (const label celli, const Type &value, const bool forceReference=false)
 Set reference level for solution. More...
 
void setReferences (const labelUList &cellLabels, const Type &value, const bool forceReference=false)
 Set reference level for solution. More...
 
void setReferences (const labelUList &cellLabels, const UList< Type > &values, const bool forceReference=false)
 Set reference level for solution. More...
 
void setComponentReference (const label patchi, const label facei, const direction cmpt, const scalar value)
 
void relax (const scalar alpha)
 Relax matrix (for steady-state solution). More...
 
void relax ()
 Relax matrix (for steady-state solution). More...
 
void boundaryManipulate (typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
 Manipulate based on a boundary field. More...
 
autoPtr< fvSolversolver (const dictionary &)
 Construct and return the solver. More...
 
autoPtr< fvSolversolver ()
 Construct and return the solver. More...
 
SolverPerformance< Type > solveSegregatedOrCoupled (const dictionary &)
 Solve segregated or coupled returning the solution statistics. More...
 
SolverPerformance< Type > solveSegregated (const dictionary &)
 Solve segregated returning the solution statistics. More...
 
SolverPerformance< Type > solveCoupled (const dictionary &)
 Solve coupled returning the solution statistics. More...
 
SolverPerformance< Type > solve (const dictionary &)
 Solve returning the solution statistics. More...
 
SolverPerformance< Type > solve ()
 Solve returning the solution statistics. More...
 
tmp< Field< Type > > residual () const
 Return the matrix residual. More...
 
tmp< scalarFieldD () const
 Return the matrix scalar diagonal. More...
 
tmp< Field< Type > > DD () const
 Return the matrix Type diagonal. More...
 
tmp< volScalarFieldA () const
 Return the central coefficient. More...
 
tmp< GeometricField< Type, fvPatchField, volMesh > > H () const
 Return the H operation source. More...
 
tmp< volScalarFieldH1 () const
 Return H(1) More...
 
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux () const
 Return the face-flux field from the matrix. More...
 
const dictionarysolverDict () const
 Return the solver dictionary taking into account finalIteration. More...
 
void operator= (const fvMatrix< Type > &)
 
void operator= (const tmp< fvMatrix< Type >> &)
 
void 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+= (const zero &)
 
void operator-= (const zero &)
 
void operator*= (const volScalarField::Internal &)
 
void operator*= (const tmp< volScalarField::Internal > &)
 
void operator*= (const tmp< volScalarField > &)
 
void operator*= (const dimensioned< scalar > &)
 
template<>
void setComponentReference (const label patchi, const label facei, const direction, const scalar value)
 
template<>
Foam::autoPtr< Foam::fvMatrix< Foam::scalar >::fvSolversolver (const dictionary &solverControls)
 
template<>
Foam::solverPerformance solveSegregated (const dictionary &solverControls)
 
template<>
Foam::tmp< Foam::scalarFieldresidual () const
 
template<>
Foam::tmp< Foam::volScalarFieldH () const
 
template<>
Foam::tmp< Foam::volScalarFieldH1 () const
 
template<>
void setComponentReference (const label patchi, const label facei, const direction, const scalar value)
 
template<>
autoPtr< fvMatrix< scalar >::fvSolversolver (const dictionary &)
 
template<>
solverPerformance solveSegregated (const dictionary &)
 
template<>
tmp< scalarFieldresidual () const
 
template<>
tmp< volScalarFieldH () const
 
template<>
tmp< volScalarFieldH1 () 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. More...
 
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. More...
 
template<class Type2 >
void subtractFromInternalField (const labelUList &addr, const tmp< Field< Type2 >> &tpf, Field< Type2 > &intf) const
 
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. More...
 

Friends

class fvSolver
 Declare friendship with the fvSolver class. More...
 
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 >> &)
 
Ostreamoperator (Ostream &, const fvMatrix< Type > &)
 

Detailed Description

template<class Type>
class Foam::fvMatrix< Type >

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.

Source files

Definition at line 76 of file fvPatchField.H.

Member Typedef Documentation

◆ surfaceTypeFieldPtr

Declare return type of the faceFluxCorrectionPtr() function.

Definition at line 338 of file fvMatrix.H.

Constructor & Destructor Documentation

◆ fvMatrix() [1/4]

fvMatrix ( const GeometricField< Type, fvPatchField, volMesh > &  psi,
const dimensionSet ds 
)

Construct given a field to solve for.

Definition at line 268 of file fvMatrix.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), DebugInFunction, Foam::endl(), forAll, psi, and Foam::Zero.

Here is the call graph for this function:

◆ fvMatrix() [2/4]

fvMatrix ( const fvMatrix< Type > &  fvm)

Copy construct.

Definition at line 319 of file fvMatrix.C.

References DebugInFunction, and Foam::endl().

Here is the call graph for this function:

◆ fvMatrix() [3/4]

fvMatrix ( const tmp< fvMatrix< Type >> &  tfvm)

Copy/move construct from tmp<fvMatrix<Type>>

Definition at line 345 of file fvMatrix.C.

References DebugInFunction, and Foam::endl().

Here is the call graph for this function:

◆ fvMatrix() [4/4]

fvMatrix ( const GeometricField< Type, fvPatchField, volMesh > &  psi,
Istream is 
)

Construct from Istream given field to solve for.

Definition at line 397 of file fvMatrix.C.

References DebugInFunction, Foam::endl(), forAll, psi, and Foam::Zero.

Here is the call graph for this function:

◆ ~fvMatrix()

~fvMatrix ( )
virtual

Destructor.

Definition at line 452 of file fvMatrix.C.

References DebugInFunction, Foam::deleteDemandDrivenData(), and Foam::endl().

Here is the call graph for this function:

Member Function Documentation

◆ addToInternalField() [1/2]

void addToInternalField ( const labelUList addr,
const Field< Type2 > &  pf,
Field< Type2 > &  intf 
) const
protected

Add patch contribution to internal field.

Definition at line 43 of file fvMatrix.C.

References Foam::abort(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, forAll, and UList< T >::size().

Here is the call graph for this function:

◆ addToInternalField() [2/2]

void addToInternalField ( const labelUList addr,
const tmp< Field< Type2 >> &  tpf,
Field< Type2 > &  intf 
) const
protected

Definition at line 67 of file fvMatrix.C.

◆ subtractFromInternalField() [1/2]

void subtractFromInternalField ( const labelUList addr,
const Field< Type2 > &  pf,
Field< Type2 > &  intf 
) const
protected

Subtract patch contribution from internal field.

Definition at line 81 of file fvMatrix.C.

References Foam::abort(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, forAll, and UList< T >::size().

Here is the call graph for this function:

◆ subtractFromInternalField() [2/2]

void subtractFromInternalField ( const labelUList addr,
const tmp< Field< Type2 >> &  tpf,
Field< Type2 > &  intf 
) const
protected

Definition at line 105 of file fvMatrix.C.

◆ addBoundaryDiag()

void addBoundaryDiag ( scalarField diag,
const direction  cmpt 
) const
protected

Definition at line 118 of file fvMatrix.C.

References Foam::component(), Foam::diag(), and forAll.

Referenced by fvMatrix< Type >::residual().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ addCmptAvBoundaryDiag()

void addCmptAvBoundaryDiag ( scalarField diag) const
protected

Definition at line 136 of file fvMatrix.C.

References Foam::cmptAv(), Foam::diag(), and forAll.

Here is the call graph for this function:

◆ addBoundarySource()

void addBoundarySource ( Field< Type > &  source,
const bool  couples = true 
) const
protected

Definition at line 152 of file fvMatrix.C.

References Foam::cmptMultiply(), fvPatchField< Type >::coupled(), forAll, and fvPatchField< Type >::patchNeighbourField().

Referenced by fvMatrix< Type >::residual().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ setValuesFromList()

void setValuesFromList ( const labelUList cellLabels,
const ListType< Type > &  values 
)
protected

Set solution in given cells to the specified values.

Definition at line 185 of file fvMatrix.C.

References cells, Foam::diag(), forAll, Foam::stringOps::lower(), mesh, primitiveFieldRef(), psi, Foam::stringOps::upper(), Foam::HashTableOps::values(), and Foam::Zero.

Here is the call graph for this function:

◆ ClassName()

ClassName ( "fvMatrix< Type >"  )

◆ clone()

Foam::tmp< Foam::fvMatrix< Type > > clone ( ) const

Clone.

Definition at line 440 of file fvMatrix.C.

◆ psi()

◆ dimensions()

const dimensionSet& dimensions ( ) const
inline

Definition at line 292 of file fvMatrix.H.

Referenced by explicitPorositySource::addSup(), interRegionExplicitPorositySource::addSup(), meanVelocityForce::addSup(), SemiImplicitSource< Type >::addSup(), directionalPressureGradientExplicitSource::addSup(), rotorDiskSource::addSup(), Foam::checkMethod(), and fvMatrix< Type >::operator*=().

Here is the caller graph for this function:

◆ source() [1/2]

◆ source() [2/2]

const Field<Type>& source ( ) const
inline

Definition at line 302 of file fvMatrix.H.

◆ internalCoeffs() [1/2]

const FieldField<Field, Type>& internalCoeffs ( ) const
inline

fvBoundary scalar field containing pseudo-matrix coeffs for internal cells

Definition at line 309 of file fvMatrix.H.

Referenced by dynamicOversetFvMesh::addInterpolation(), gaussConvectionScheme< Type >::fvmDiv(), gaussLaplacianScheme< Type, GType >::fvmLaplacianUncorrected(), dynamicOversetFvMesh::normalisation(), dynamicOversetFvMesh::solve(), and dynamicOversetFvMesh::write().

Here is the caller graph for this function:

◆ internalCoeffs() [2/2]

FieldField<Field, Type>& internalCoeffs ( )
inline

fvBoundary scalar field containing pseudo-matrix coeffs for internal cells

Definition at line 316 of file fvMatrix.H.

◆ boundaryCoeffs() [1/2]

const FieldField<Field, Type>& boundaryCoeffs ( ) const
inline

fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells

Definition at line 323 of file fvMatrix.H.

Referenced by dynamicOversetFvMesh::addInterpolation(), gaussConvectionScheme< Type >::fvmDiv(), gaussLaplacianScheme< Type, GType >::fvmLaplacianUncorrected(), dynamicOversetFvMesh::solve(), and dynamicOversetFvMesh::write().

Here is the caller graph for this function:

◆ boundaryCoeffs() [2/2]

FieldField<Field, Type>& boundaryCoeffs ( )
inline

fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells

Definition at line 330 of file fvMatrix.H.

◆ faceFluxCorrectionPtr()

surfaceTypeFieldPtr& faceFluxCorrectionPtr ( )
inline

Return pointer to face-flux non-orthogonal correction field.

Definition at line 341 of file fvMatrix.H.

Referenced by gaussLaplacianScheme< Type, GType >::fvmLaplacian().

Here is the caller graph for this function:

◆ setValues() [1/3]

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 465 of file fvMatrix.C.

Referenced by FixedValueConstraint< Type >::constrain(), fixedTemperatureConstraint::constrain(), fixedInternalValueFvPatchField< Type >::manipulateMatrix(), epsilonWallFunctionFvPatchScalarField::manipulateMatrix(), and omegaWallFunctionFvPatchScalarField::manipulateMatrix().

Here is the caller graph for this function:

◆ setValues() [2/3]

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 476 of file fvMatrix.C.

References Foam::HashTableOps::values().

Here is the call graph for this function:

◆ setValues() [3/3]

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 487 of file fvMatrix.C.

References Foam::HashTableOps::values().

Here is the call graph for this function:

◆ setReference()

void setReference ( const label  celli,
const Type &  value,
const bool  forceReference = false 
)

Set reference level for solution.

Definition at line 498 of file fvMatrix.C.

References Foam::diag().

Referenced by Foam::CorrectPhi(), simple::mainIter(), and adjointSimple::mainIter().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ setReferences() [1/2]

void setReferences ( const labelUList cellLabels,
const Type &  value,
const bool  forceReference = false 
)

Set reference level for solution.

Definition at line 514 of file fvMatrix.C.

References cellId, Foam::diag(), and forAll.

Here is the call graph for this function:

◆ setReferences() [2/2]

void setReferences ( const labelUList cellLabels,
const UList< Type > &  values,
const bool  forceReference = false 
)

Set reference level for solution.

Definition at line 537 of file fvMatrix.C.

References cellId, Foam::diag(), forAll, and Foam::HashTableOps::values().

Here is the call graph for this function:

◆ setComponentReference() [1/3]

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 38 of file fvMatrixSolve.C.

References Foam::diag().

Here is the call graph for this function:

◆ relax() [1/2]

void relax ( const scalar  alpha)

◆ relax() [2/2]

void relax ( )

Relax matrix (for steady-state solution).

alpha is read from controlDict

Definition at line 707 of file fvMatrix.C.

References Foam::name(), and relax().

Here is the call graph for this function:

◆ boundaryManipulate()

void boundaryManipulate ( typename GeometricField< Type, fvPatchField, volMesh >::Boundary &  values)

Manipulate based on a boundary field.

Definition at line 724 of file fvMatrix.C.

References forAll.

Referenced by kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), adjointSimple::mainIter(), and adjointMeshMovementSolver::solve().

Here is the caller graph for this function:

◆ solver() [1/4]

Construct and return the solver.

Use the given solver controls

◆ solver() [2/4]

Foam::autoPtr< typename Foam::fvMatrix< Type >::fvSolver > solver ( )

Construct and return the solver.

Solver controls read from fvSolution

Definition at line 311 of file fvMatrixSolve.C.

◆ solveSegregatedOrCoupled()

Foam::SolverPerformance< Type > solveSegregatedOrCoupled ( const dictionary solverControls)

Solve segregated or coupled returning the solution statistics.

Use the given solver controls

Definition at line 62 of file fvMatrixSolve.C.

References addProfiling, Foam::expressions::patchExpr::debug, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInFunction, dictionary::getOrDefault(), Foam::Info, messageStream::masterStream(), mesh, dictionary::readIfPresent(), regionName, and solve().

Referenced by velocityComponentLaplacianFvMotionSolver::solve(), velocityLaplacianFvMotionSolver::solve(), displacementComponentLaplacianFvMotionSolver::solve(), displacementSBRStressFvMotionSolver::solve(), solidBodyDisplacementLaplacianFvMotionSolver::solve(), displacementLaplacianFvMotionSolver::solve(), and fvMesh::solve().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ solveSegregated() [1/3]

Foam::SolverPerformance< Type > solveSegregated ( const dictionary solverControls)

Solve segregated returning the solution statistics.

Use the given solver controls

Definition at line 115 of file fvMatrixSolve.C.

References Field< Type >::component(), Foam::expressions::patchExpr::debug, Foam::diag(), Foam::endl(), Foam::Info, messageStream::masterStream(), mesh, Foam::New(), SolverPerformance< Type >::print(), psi, refPtr< Container< Type > >::ref(), solve(), and SolverPerformance< Type >::solverName().

Here is the call graph for this function:

◆ solveCoupled()

Foam::SolverPerformance< Type > solveCoupled ( const dictionary solverControls)

Solve coupled returning the solution statistics.

Use the given solver controls

Definition at line 241 of file fvMatrixSolve.C.

References Foam::expressions::patchExpr::debug, Foam::diag(), LduMatrix< Type, DType, LUType >::diag(), Foam::endl(), Foam::Info, Foam::stringOps::lower(), messageStream::masterStream(), mesh, SolverPerformance< Type >::print(), psi, and Foam::stringOps::upper().

Here is the call graph for this function:

◆ solve() [1/2]

◆ solve() [2/2]

Foam::SolverPerformance< Type > solve ( )

Solve returning the solution statistics.

Solver controls read from fvSolution

Definition at line 325 of file fvMatrixSolve.C.

References fvMatrix< Type >::solverDict().

Here is the call graph for this function:

◆ residual() [1/3]

Foam::tmp< Foam::Field< Type > > residual ( ) const

Return the matrix residual.

Definition at line 332 of file fvMatrixSolve.C.

References fvMatrix< Type >::addBoundaryDiag(), fvMatrix< Type >::addBoundarySource(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), Field< Type >::component(), GeometricField< Type, PatchField, GeoMesh >::primitiveField(), tmp< T >::ref(), Field< Type >::replace(), lduMatrix::residual(), and Foam::Zero.

Here is the call graph for this function:

◆ D()

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

Return the matrix scalar diagonal.

Definition at line 737 of file fvMatrix.C.

References Foam::diag(), and tmp< T >::ref().

Here is the call graph for this function:

◆ DD()

Foam::tmp< Foam::Field< Type > > DD ( ) const

Return the matrix Type diagonal.

Definition at line 746 of file fvMatrix.C.

References fvPatchField< Type >::coupled(), Foam::diag(), forAll, and tmp< T >::ref().

Here is the call graph for this function:

◆ A()

Return the central coefficient.

Definition at line 770 of file fvMatrix.C.

References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), D, Foam::dimVol, IOobject::NO_READ, IOobject::NO_WRITE, GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), and tmp< T >::ref().

Referenced by meanVelocityForce::constrain(), directionalPressureGradientExplicitSource::constrain(), incompressiblePrimalSolver::correctBoundaryConditions(), and simple::mainIter().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ H() [1/3]

Return the H operation source.

Definition at line 799 of file fvMatrix.C.

References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::dimVol, lduMatrix::H(), IOobject::NO_READ, IOobject::NO_WRITE, GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), tmp< T >::ref(), Field< Type >::replace(), GeometricField< Type, PatchField, GeoMesh >::replace(), and Foam::Zero.

Referenced by incompressiblePrimalSolver::correctBoundaryConditions(), and simple::mainIter().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ H1() [1/3]

Return H(1)

Definition at line 861 of file fvMatrix.C.

References Foam::component(), fvPatchField< Type >::coupled(), Foam::dimVol, forAll, lduMatrix::H1(), IOobject::NO_READ, IOobject::NO_WRITE, and GeometricField< Type, PatchField, GeoMesh >::ref().

Referenced by simple::mainIter().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ flux()

Return the face-flux field from the matrix.

Definition at line 909 of file fvMatrix.C.

References Foam::abort(), GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), Foam::cmptMultiply(), lduMatrix::faceH(), Foam::FatalError, FatalErrorInFunction, forAll, IOobject::NO_READ, IOobject::NO_WRITE, GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), and tmp< T >::ref().

Referenced by Implicit< CloudType >::cacheFields(), incompressiblePrimalSolver::correctBoundaryConditions(), Foam::CorrectPhi(), scalarTransport::execute(), simple::mainIter(), adjointSimple::mainIter(), twoPhaseSystem::solve(), kinematicSingleLayer::solveThickness(), and MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::solveYi().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ solverDict()

const Foam::dictionary & solverDict ( ) const

Return the solver dictionary taking into account finalIteration.

Definition at line 996 of file fvMatrix.C.

Referenced by velocityComponentLaplacianFvMotionSolver::solve(), velocityLaplacianFvMotionSolver::solve(), displacementComponentLaplacianFvMotionSolver::solve(), displacementSBRStressFvMotionSolver::solve(), solidBodyDisplacementLaplacianFvMotionSolver::solve(), displacementLaplacianFvMotionSolver::solve(), and fvMatrix< Type >::solve().

Here is the caller graph for this function:

◆ operator=() [1/2]

void operator= ( const fvMatrix< Type > &  fvmv)

Definition at line 1012 of file fvMatrix.C.

References Foam::abort(), Foam::FatalError, FatalErrorInFunction, and lduMatrix::operator=().

Here is the call graph for this function:

◆ operator=() [2/2]

void operator= ( const tmp< fvMatrix< Type >> &  tfvmv)

Definition at line 1046 of file fvMatrix.C.

◆ negate()

void negate ( )

Definition at line 1054 of file fvMatrix.C.

References lduMatrix::negate().

Here is the call graph for this function:

◆ operator+=() [1/7]

void operator+= ( const fvMatrix< Type > &  fvmv)

Definition at line 1069 of file fvMatrix.C.

References Foam::checkMethod(), and lduMatrix::operator+=().

Here is the call graph for this function:

◆ operator+=() [2/7]

void operator+= ( const tmp< fvMatrix< Type >> &  tfvmv)

Definition at line 1095 of file fvMatrix.C.

◆ operator-=() [1/7]

void operator-= ( const fvMatrix< Type > &  fvmv)

Definition at line 1103 of file fvMatrix.C.

References Foam::checkMethod(), and lduMatrix::operator-=().

Here is the call graph for this function:

◆ operator-=() [2/7]

void operator-= ( const tmp< fvMatrix< Type >> &  tfvmv)

Definition at line 1127 of file fvMatrix.C.

◆ operator+=() [3/7]

void operator+= ( const DimensionedField< Type, volMesh > &  su)

Definition at line 1136 of file fvMatrix.C.

References Foam::checkMethod().

Here is the call graph for this function:

◆ operator+=() [4/7]

void operator+= ( const tmp< DimensionedField< Type, volMesh >> &  tsu)

Definition at line 1147 of file fvMatrix.C.

◆ operator+=() [5/7]

void operator+= ( const tmp< GeometricField< Type, fvPatchField, volMesh >> &  tsu)

Definition at line 1158 of file fvMatrix.C.

◆ operator-=() [3/7]

void operator-= ( const DimensionedField< Type, volMesh > &  su)

Definition at line 1169 of file fvMatrix.C.

References Foam::checkMethod().

Here is the call graph for this function:

◆ operator-=() [4/7]

void operator-= ( const tmp< DimensionedField< Type, volMesh >> &  tsu)

Definition at line 1180 of file fvMatrix.C.

◆ operator-=() [5/7]

void operator-= ( const tmp< GeometricField< Type, fvPatchField, volMesh >> &  tsu)

Definition at line 1191 of file fvMatrix.C.

◆ operator+=() [6/7]

void operator+= ( const dimensioned< Type > &  su)

Definition at line 1202 of file fvMatrix.C.

References psi.

◆ operator-=() [6/7]

void operator-= ( const dimensioned< Type > &  su)

Definition at line 1212 of file fvMatrix.C.

References psi.

◆ operator+=() [7/7]

void operator+= ( const zero )

Definition at line 1222 of file fvMatrix.C.

◆ operator-=() [7/7]

void operator-= ( const zero )

Definition at line 1230 of file fvMatrix.C.

◆ operator*=() [1/4]

void operator*= ( const volScalarField::Internal dsf)

Definition at line 1238 of file fvMatrix.C.

References Foam::abort(), fvMatrix< Type >::dimensions(), Foam::FatalError, FatalErrorInFunction, forAll, and lduMatrix::operator*=().

Here is the call graph for this function:

◆ operator*=() [2/4]

void operator*= ( const tmp< volScalarField::Internal > &  tdsf)

Definition at line 1268 of file fvMatrix.C.

◆ operator*=() [3/4]

void operator*= ( const tmp< volScalarField > &  tvsf)

Definition at line 1279 of file fvMatrix.C.

◆ operator*=() [4/4]

void operator*= ( const dimensioned< scalar > &  ds)

Definition at line 1290 of file fvMatrix.C.

References lduMatrix::operator*=().

Here is the call graph for this function:

◆ setComponentReference() [2/3]

void setComponentReference ( const label  patchi,
const label  facei,
const  direction,
const scalar  value 
)

Definition at line 38 of file fvScalarMatrix.C.

References Foam::diag().

Here is the call graph for this function:

◆ solver() [3/4]

Foam::autoPtr< Foam::fvMatrix< Foam::scalar >::fvSolver > solver ( const dictionary solverControls)

Definition at line 63 of file fvScalarMatrix.C.

References addProfiling, Foam::expressions::patchExpr::debug, Foam::diag(), Foam::endl(), Foam::Info, messageStream::masterStream(), mesh, Foam::New(), regionName, and solve().

Here is the call graph for this function:

◆ solveSegregated() [2/3]

Foam::solverPerformance solveSegregated ( const dictionary solverControls)

Definition at line 150 of file fvScalarMatrix.C.

References Foam::expressions::patchExpr::debug, Foam::diag(), Foam::endl(), Foam::Info, messageStream::masterStream(), mesh, Foam::New(), SolverPerformance< Type >::print(), psi, and solve().

Here is the call graph for this function:

◆ residual() [2/3]

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

Definition at line 199 of file fvScalarMatrix.C.

References psi, refPtr< Container< Type > >::ref(), and Foam::Zero.

Here is the call graph for this function:

◆ H() [2/3]

Definition at line 228 of file fvScalarMatrix.C.

References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::dimVol, H(), GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), and tmp< T >::ref().

Here is the call graph for this function:

◆ H1() [2/3]

Definition at line 260 of file fvScalarMatrix.C.

References Foam::dimVol, and GeometricField< Type, PatchField, GeoMesh >::ref().

Here is the call graph for this function:

◆ setComponentReference() [3/3]

void setComponentReference ( const label  patchi,
const label  facei,
const  direction,
const scalar  value 
)

◆ solver() [4/4]

autoPtr< fvMatrix< scalar >::fvSolver > solver ( const dictionary )

◆ solveSegregated() [3/3]

solverPerformance solveSegregated ( const dictionary )

◆ residual() [3/3]

tmp< scalarField > residual ( ) const

◆ H() [3/3]

tmp< volScalarField > H ( ) const

◆ H1() [3/3]

tmp< volScalarField > H1 ( ) const

Friends And Related Function Documentation

◆ fvSolver

friend class fvSolver
friend

Declare friendship with the fvSolver class.

Definition at line 148 of file fvMatrix.H.

◆ operator& [1/4]

tmp<GeometricField<Type, fvPatchField, volMesh> > operator& ( const fvMatrix< Type > &  ,
const DimensionedField< Type, volMesh > &   
)
friend

◆ operator& [2/4]

tmp<GeometricField<Type, fvPatchField, volMesh> > operator& ( const fvMatrix< Type > &  ,
const tmp< GeometricField< Type, fvPatchField, volMesh >> &   
)
friend

◆ operator& [3/4]

tmp<GeometricField<Type, fvPatchField, volMesh> > operator& ( const tmp< fvMatrix< Type >> &  ,
const DimensionedField< Type, volMesh > &   
)
friend

◆ operator& [4/4]

tmp<GeometricField<Type, fvPatchField, volMesh> > operator& ( const tmp< fvMatrix< Type >> &  ,
const tmp< GeometricField< Type, fvPatchField, volMesh >> &   
)
friend

◆ operator

Ostream& operator ( Ostream ,
const fvMatrix< Type > &   
)
friend

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