43 if (!isA<wallFvPatch>(patch()))
46 <<
"Invalid wall function specification" <<
nl
47 <<
" Patch type for patch " << patch().name()
48 <<
" must be wall" <<
nl
49 <<
" Current patch type is " << patch().type() <<
nl <<
endl
75 fixedValueFvPatchScalarField(
p, iF,
dict),
76 Cmu_(
dict.lookupOrDefault<scalar>(
"Cmu",
dimless, 0.09)),
77 kappa_(
dict.lookupOrDefault<scalar>(
"kappa",
dimless, 0.41)),
78 E_(
dict.lookupOrDefault<scalar>(
"E",
dimless, 9.8)),
79 yPlusLam_(yPlusLam(kappa_, E_))
93 fixedValueFvPatchScalarField(ptf,
p, iF, mapper, false),
97 yPlusLam_(ptf.yPlusLam_)
100 mapper(*
this, ptf, [&](){
return this->patchInternalField(); });
110 fixedValueFvPatchScalarField(wfpsf, iF),
112 kappa_(wfpsf.kappa_),
114 yPlusLam_(wfpsf.yPlusLam_)
130 refCast<const nutWallFunctionFvPatchScalarField>
132 turbModel.
nut()().boundaryField()[
patchi]
143 const scalar ypl0 = 12.0;
147 for (
int i=0; i<10; i++)
159 for (
int i=0; i<10; i++)
181 mapper(*
this, ptf, [&](){
return this->patchInternalField(); });
194 fixedValueFvPatchScalarField::updateCoeffs();
201 writeLocalEntries(os);
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Abstract base class for field mapping.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
virtual void write(Ostream &) const
Write.
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
scalar yPlusLam() const
Return the Y+ at the edge of the laminar sublayer.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
virtual void checkType()
Check the type of the patch.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
errorManip< error > abort(error &err)
const dimensionSet dimless
dimensionedScalar log(const dimensionedScalar &ds)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)