Go to the documentation of this file.
44 if (!isA<wallFvPatch>(
patch()))
47 <<
"Invalid wall function specification" <<
nl
48 <<
" Patch type for patch " <<
patch().
name()
49 <<
" must be wall" <<
nl
50 <<
" Current patch type is " <<
patch().type() <<
nl <<
endl
72 for (
int i=0; i<10; i++)
93 yPlusLam_(yPlusLam(kappa_, E_))
125 Cmu_(
dict.lookupOrDefault<scalar>(
"Cmu", 0.09)),
126 kappa_(
dict.lookupOrDefault<scalar>(
"kappa", 0.41)),
127 E_(
dict.lookupOrDefault<scalar>(
"E", 9.8)),
128 yPlusLam_(yPlusLam(kappa_, E_))
184 const v2fBase& v2fModel = refCast<const v2fBase>(turbModel);
216 scalar v2c = v2[faceCellI];
217 scalar epsc =
epsilon[faceCellI];
218 scalar kc =
k[faceCellI];
220 f[faceI] =
N*v2c*epsc/(
sqr(kc) + ROOTVSMALL);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
This boundary condition provides a turbulence damping function, f, wall function condition for low- a...
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
static word groupName(Name name, const word &group)
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field, sets Updated to false.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
const DimensionedField< Type, volMesh > & dimensionedInternalField() const
Return dimensioned internal field reference.
static const word propertiesName
Default name of the turbulence properties dictionary.
scalar Cmu_
Cmu coefficient.
const char *const group
Group name for atomic constants.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const nearWallDist & y() const
Return the near wall distances.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
bool updated() const
Return true if the boundary condition has already been updated.
virtual void write(Ostream &) const
Write.
dimensionedScalar pow025(const dimensionedScalar &ds)
const word & name() const
Return name.
fWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual void checkType()
Check the type of the patch.
Abstract base class for turbulence models (RAS, LES and laminar).
A list of keyword definitions, which are a keyword followed by any number of values (e....
label index() const
Return the index of this patch in the fvBoundaryMesh.
dimensionedScalar log(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
scalar kappa_
Von Karman constant.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
errorManip< error > abort(error &err)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
commsTypes
Types of communications.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
makePatchTypeField(fvPatchScalarField, fWallFunctionFvPatchScalarField)
const objectRegistry & db() const
Return local objectRegistry.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
This function object evaluates and outputs turbulence y+ for turbulence models. The field is stored o...
scalar yPlusLam(const scalar kappa, const scalar E)
Calculate the Y+ at the edge of the laminar sublayer.
label k
Boltzmann constant.
virtual void evaluate(const Pstream::commsTypes)
Evaluate the patchField.
virtual const labelUList & faceCells() const
Return faceCells.
Foam::fvPatchFieldMapper.
const fvPatch & patch() const
Return patch.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Generic GeometricField class.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Abstract base-class for v2-f models to provide BCs access to the v2 and f fields.
virtual void write(Ostream &) const
Write.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.