fWallFunctionFvPatchScalarField Class Reference

This boundary condition provides a wall constraint on the elliptic relaxation factor, f, which is executed in the v2-f based turbulence closure models for low- and high-Reynolds number turbulent flow cases. More...

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

Public Member Functions

 TypeName ("fWallFunction")
 Runtime type information. More...
 
 fWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Construct from patch and internal field. More...
 
 fWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 Construct from patch, internal field and dictionary. More...
 
 fWallFunctionFvPatchScalarField (const fWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 
 fWallFunctionFvPatchScalarField (const fWallFunctionFvPatchScalarField &)
 Construct as copy. More...
 
virtual tmp< fvPatchScalarFieldclone () const
 Construct and return a clone. More...
 
 fWallFunctionFvPatchScalarField (const fWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &)
 Construct as copy setting internal field reference. More...
 
virtual tmp< fvPatchScalarFieldclone (const DimensionedField< scalar, volMesh > &iF) const
 Construct and return a clone setting internal field reference. More...
 
virtual void updateCoeffs ()
 Update the coefficients associated with the patch field. More...
 
- Public Member Functions inherited from fixedValueFvPatchField< scalar >
 TypeName ("fixedValue")
 Runtime type information. More...
 
 fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Construct from patch and internal field. More...
 
 fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const scalar &value)
 Construct from patch, internal field and value. More...
 
 fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &, const bool valueRequired=true)
 Construct from patch, internal field and dictionary. More...
 
 fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 Construct by mapping the given fixedValueFvPatchField<Type> More...
 
 fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &)
 Construct as copy. More...
 
 fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &, const DimensionedField< scalar, volMesh > &)
 Construct as copy setting internal field reference. More...
 
virtual bool fixesValue () const
 Return true if this patch field fixes a value. More...
 
virtual bool assignable () const
 Return false: this patch field is not altered by assignment. More...
 
virtual tmp< Field< scalar > > valueInternalCoeffs (const tmp< scalarField > &) const
 Return the matrix diagonal coefficients corresponding to the. More...
 
virtual tmp< Field< scalar > > valueBoundaryCoeffs (const tmp< scalarField > &) const
 Return the matrix source coefficients corresponding to the. More...
 
virtual tmp< Field< scalar > > gradientInternalCoeffs () const
 Return the matrix diagonal coefficients corresponding to the. More...
 
virtual tmp< Field< scalar > > gradientBoundaryCoeffs () const
 Return the matrix source coefficients corresponding to the. More...
 
virtual void write (Ostream &) const
 Write. More...
 
virtual void operator= (const UList< scalar > &)
 
virtual void operator= (const fvPatchField< scalar > &)
 
virtual void operator= (const scalar &)
 
virtual void operator+= (const fvPatchField< scalar > &)
 
virtual void operator+= (const Field< scalar > &)
 
virtual void operator+= (const scalar &)
 
virtual void operator-= (const fvPatchField< scalar > &)
 
virtual void operator-= (const Field< scalar > &)
 
virtual void operator-= (const scalar &)
 
virtual void operator*= (const fvPatchField< scalar > &)
 
virtual void operator*= (const Field< scalar > &)
 
virtual void operator*= (const scalar)
 
virtual void operator/= (const fvPatchField< scalar > &)
 
virtual void operator/= (const Field< scalar > &)
 
virtual void operator/= (const scalar)
 
- Public Member Functions inherited from fvPatchField< scalar >
 TypeName ("fvPatchField")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (tmp, fvPatchField, patch,(const fvPatch &p, const DimensionedField< scalar, volMesh > &iF),(p, iF))
 
 declareRunTimeSelectionTable (tmp, fvPatchField, patchMapper,(const fvPatchField< scalar > &ptf, const fvPatch &p, const DimensionedField< scalar, volMesh > &iF, const fvPatchFieldMapper &m),(dynamic_cast< const fvPatchFieldType & >(ptf), p, iF, m))
 
 declareRunTimeSelectionTable (tmp, fvPatchField, dictionary,(const fvPatch &p, const DimensionedField< scalar, volMesh > &iF, const dictionary &dict),(p, iF, dict))
 
 fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Construct from patch and internal field. More...
 
 fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const scalar &value)
 Construct from patch, internal field and value. More...
 
 fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const word &patchType)
 Construct from patch and internal field and patch type. More...
 
 fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const Field< scalar > &)
 Construct from patch and internal field and patch field. More...
 
 fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &, const bool valueRequired=true)
 Construct from patch, internal field and dictionary. More...
 
 fvPatchField (const fvPatchField< scalar > &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 Construct by mapping the given fvPatchField onto a new patch. More...
 
 fvPatchField (const fvPatchField< scalar > &)
 Construct as copy. More...
 
 fvPatchField (const fvPatchField< scalar > &, const DimensionedField< scalar, volMesh > &)
 Construct as copy setting internal field reference. More...
 
Foam::tmp< Foam::fvPatchField< scalar > > NewCalculatedType (const fvPatchField< Type2 > &pf)
 
virtual ~fvPatchField ()=default
 Destructor. More...
 
virtual bool coupled () const
 Return true if this patch field is coupled. More...
 
const objectRegistrydb () const
 Return local objectRegistry. More...
 
const fvPatchpatch () const
 Return patch. More...
 
const DimensionedField< scalar, volMesh > & internalField () const
 Return dimensioned internal field reference. More...
 
const Field< scalar > & primitiveField () const
 Return internal field reference. More...
 
const wordpatchType () const
 Optional patch type. More...
 
wordpatchType ()
 Optional patch type. More...
 
bool updated () const
 Return true if the boundary condition has already been updated. More...
 
bool manipulatedMatrix () const
 Return true if the matrix has already been manipulated. More...
 
virtual void autoMap (const fvPatchFieldMapper &)
 Map (and resize as needed) from self given a mapping object. More...
 
virtual void rmap (const fvPatchField< scalar > &, const labelList &)
 Reverse map the given fvPatchField onto this fvPatchField. More...
 
virtual tmp< Field< scalar > > snGrad () const
 Return patch-normal gradient. More...
 
virtual tmp< Field< scalar > > snGrad (const scalarField &deltaCoeffs) const
 Return patch-normal gradient for coupled-patches. More...
 
virtual void updateWeightedCoeffs (const scalarField &weights)
 Update the coefficients associated with the patch field. More...
 
virtual tmp< Field< scalar > > patchInternalField () const
 Return internal field next to patch as patch field. More...
 
virtual void patchInternalField (Field< scalar > &) const
 Return internal field next to patch as patch field. More...
 
virtual tmp< Field< scalar > > patchNeighbourField () const
 Return patchField on the opposite patch of a coupled patch. More...
 
virtual void initEvaluate (const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
 Initialise the evaluation of the patch field. More...
 
virtual void evaluate (const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
 Evaluate the patch field, sets Updated to false. More...
 
virtual tmp< Field< scalar > > valueInternalCoeffs (const tmp< Field< scalar >> &) const
 Return the matrix diagonal coefficients corresponding to the. More...
 
virtual tmp< Field< scalar > > valueBoundaryCoeffs (const tmp< Field< scalar >> &) const
 Return the matrix source coefficients corresponding to the. More...
 
virtual tmp< Field< scalar > > gradientInternalCoeffs (const scalarField &deltaCoeffs) const
 Return the matrix diagonal coefficients corresponding to the. More...
 
virtual tmp< Field< scalar > > gradientBoundaryCoeffs (const scalarField &deltaCoeffs) const
 Return the matrix source coefficients corresponding to the. More...
 
virtual void manipulateMatrix (fvMatrix< scalar > &matrix)
 Manipulate matrix. More...
 
virtual void manipulateMatrix (fvMatrix< scalar > &matrix, const scalarField &weights)
 Manipulate matrix with given weights. More...
 
void check (const fvPatchField< scalar > &) const
 Check fvPatchField<Type> against given fvPatchField<Type> More...
 
virtual void operator== (const fvPatchField< scalar > &)
 
virtual void operator== (const Field< scalar > &)
 
virtual void operator== (const scalar &)
 

Additional Inherited Members

- Public Types inherited from fvPatchField< scalar >
typedef fvPatch Patch
 
typedef calculatedFvPatchField< scalar > Calculated
 
- Static Public Member Functions inherited from fvPatchField< scalar >
static tmp< fvPatchField< scalar > > New (const word &, const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Return a pointer to a new patchField created on freestore given. More...
 
static tmp< fvPatchField< scalar > > New (const word &, const word &actualPatchType, const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Return a pointer to a new patchField created on freestore given. More...
 
static tmp< fvPatchField< scalar > > New (const fvPatchField< scalar > &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 Return a pointer to a new patchField created on freestore from. More...
 
static tmp< fvPatchField< scalar > > New (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 Return a pointer to a new patchField created on freestore. More...
 
static tmp< fvPatchField< scalar > > NewCalculatedType (const fvPatch &)
 Return a pointer to a new calculatedFvPatchField created on. More...
 
static tmp< fvPatchField< scalar > > NewCalculatedType (const fvPatchField< Type2 > &)
 Return a pointer to a new calculatedFvPatchField created on. More...
 
static const wordcalculatedType ()
 Return the type of the calculated for of fvPatchField. More...
 
- Static Public Attributes inherited from fvPatchField< scalar >
static int disallowGenericFvPatchField
 Debug switch to disallow the use of genericFvPatchField. More...
 

Detailed Description

This boundary condition provides a wall constraint on the elliptic relaxation factor, f, which is executed in the v2-f based turbulence closure models for low- and high-Reynolds number turbulent flow cases.

Reference:

Remark on the blending approach for f (tag:PH):
    Popovac, M., & Hanjalić, K. (2007).
    Compound wall treatment for RANS computation of complex
    turbulent flows and heat transfer.
    Flow, turbulence and combustion, 78(2), 177-202.
    DOI:10.1007/s10494-006-9067-x

Wall-boundary expression for f in kEpsilonPhitF model (tag:LUU):
    Laurence, D. R., Uribe, J. C., & Utyuzhnikov, S. V. (2005).
    A robust formulation of the v2−f model.
    Flow, Turbulence and Combustion, 73(3-4), 169–185.
    DOI:10.1007/s10494-005-1974-8
Usage
Example of the boundary condition specification:
<patchName>
{
    // Mandatory entries
    type            fWallFunction;

    // No optional entry
}
Note
The coefficients Cmu, kappa, and E are obtained from the specified nutWallFunction in order to ensure that each patch possesses the same set of values for these coefficients.

For f, the viscous and inertial sublayer blending approaches were claimed to be inviable (PH:p. 194). Therefore, the only blending mode for the v2-f model is the stepwise mode where the viscous and inertial sublayer contributions switch over a sublayer-intersection value of y+ estimated from the kappa and E.

See also
Foam::fixedValueFvPatchField
Source files

Definition at line 101 of file fWallFunctionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ fWallFunctionFvPatchScalarField() [1/5]

fWallFunctionFvPatchScalarField ( const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF 
)

Construct from patch and internal field.

Definition at line 45 of file fWallFunctionFvPatchScalarField.C.

Referenced by fWallFunctionFvPatchScalarField::clone().

Here is the caller graph for this function:

◆ fWallFunctionFvPatchScalarField() [2/5]

fWallFunctionFvPatchScalarField ( const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF,
const dictionary dict 
)

Construct from patch, internal field and dictionary.

Definition at line 67 of file fWallFunctionFvPatchScalarField.C.

◆ fWallFunctionFvPatchScalarField() [3/5]

fWallFunctionFvPatchScalarField ( const fWallFunctionFvPatchScalarField ptf,
const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF,
const fvPatchFieldMapper mapper 
)

Construct by mapping given fWallFunctionFvPatchScalarField onto a new patch

Definition at line 55 of file fWallFunctionFvPatchScalarField.C.

◆ fWallFunctionFvPatchScalarField() [4/5]

Construct as copy.

Definition at line 78 of file fWallFunctionFvPatchScalarField.C.

◆ fWallFunctionFvPatchScalarField() [5/5]

Construct as copy setting internal field reference.

Definition at line 87 of file fWallFunctionFvPatchScalarField.C.

Member Function Documentation

◆ TypeName()

TypeName ( "fWallFunction"  )

Runtime type information.

◆ clone() [1/2]

virtual tmp<fvPatchScalarField> clone ( ) const
inlinevirtual

Construct and return a clone.

Reimplemented from fixedValueFvPatchField< scalar >.

Definition at line 145 of file fWallFunctionFvPatchScalarField.H.

References fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField().

Here is the call graph for this function:

◆ clone() [2/2]

virtual tmp<fvPatchScalarField> clone ( const DimensionedField< scalar, volMesh > &  iF) const
inlinevirtual

Construct and return a clone setting internal field reference.

Reimplemented from fixedValueFvPatchField< scalar >.

Definition at line 162 of file fWallFunctionFvPatchScalarField.H.

References fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField().

Here is the call graph for this function:

◆ updateCoeffs()


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