Go to the documentation of this file.
51 (KsPlus - 2.25)/87.75 + Cs*KsPlus,
52 sin(0.4258*(
log(KsPlus) - 0.811))
57 return (1.0 + Cs*KsPlus);
87 label faceCellI = patch().faceCells()[faceI];
89 scalar uStar = Cmu25*
sqrt(
k[faceCellI]);
90 scalar
yPlus = uStar*
y[faceI]/nuw[faceI];
91 scalar KsPlus = uStar*
Ks_[faceI]/nuw[faceI];
99 scalar limitingNutw =
max(nutw[faceI], nuw[faceI]);
117 <<
", KsPlus = " << KsPlus
118 <<
", Edash = " << Edash
119 <<
", nutw = " << nutw[faceI]
151 Ks_(ptf.
Ks_, mapper),
164 Ks_(
"Ks",
dict,
p.size()),
165 Cs_(
"Cs",
dict,
p.size())
199 nutkWallFunctionFvPatchScalarField::autoMap(m);
211 nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
214 refCast<const nutkRoughWallFunctionFvPatchScalarField>(ptf);
216 Ks_.rmap(nrwfpsf.
Ks_, addr);
217 Cs_.rmap(nrwfpsf.
Cs_, addr);
227 writeEntry(
"value", os);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
static word groupName(Name name, const word &group)
virtual void write(Ostream &) const
Write.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
dimensionedScalar sin(const dimensionedScalar &ds)
rDeltaT dimensionedInternalField()
static const word propertiesName
Default name of the turbulence properties dictionary.
const char *const group
Group name for atomic constants.
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
const nearWallDist & y() const
Return the near wall distances.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
scalar kappa_
Von Karman constant.
dimensionedScalar pow025(const dimensionedScalar &ds)
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.
nutkRoughWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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....
dimensionedScalar log(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
scalarField Ks_
Roughness height.
scalar Cmu_
Cmu coefficient.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
const double e
Elementary charge.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
virtual scalar fnRough(const scalar KsPlus, const scalar Cs) const
Compute the roughness function.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions ...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedScalar sqrt(const dimensionedScalar &ds)
This function object evaluates and outputs turbulence y+ for turbulence models. The field is stored o...
label k
Boltzmann constant.
scalarField Cs_
Roughness constant.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Generic GeometricField class.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...