Go to the documentation of this file.
69 for (
const scalar z : z0)
74 <<
"z0 field can only contain positive values. "
75 <<
"Please check input field z0."
88 const scalar w = cornerWeights[facei];
91 w*
sqrt(
k[celli])/(Cmu25*nutw.
kappa()*(
y[facei] + z0[facei]));
95 *(nutw[facei] + nuw[facei])
98 /(nutw.
kappa()*(
y[facei] + z0[facei]));
127 z0_(ptf.z0_.clone(
p.
patch()))
188 refCast<const atmOmegaWallFunctionFvPatchScalarField>(ptf);
190 z0_->rmap(atmpsf.
z0_(), addr);
211 atmOmegaWallFunctionFvPatchScalarField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
fvPatchField< scalar > fvPatchScalarField
scalar timeOutputValue() const
virtual tmp< Field< Type > > snGrad() const
virtual void write(Ostream &) const
atmOmegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
virtual tmp< volScalarField > nu() const =0
This boundary condition provides a wall constraint on the specific dissipation rate,...
A class for managing temporary objects.
This boundary condition provides a wall constraint on the specific dissipation rate (i....
const nearWallDist & y() const
virtual void rmap(const fvPatchScalarField &, const labelList &)
dimensionedScalar pow025(const dimensionedScalar &ds)
autoPtr< PatchFunction1< scalar > > z0_
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Generic templated field type.
const dimensionedScalar G0
virtual void write(Ostream &) const
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Abstract base class for turbulence models (RAS, LES and laminar).
virtual void rmap(const fvPatchField< Type > &, const labelList &)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
virtual tmp< volScalarField > k() const =0
virtual void autoMap(const fvPatchFieldMapper &)
errorManipArg< error, int > exit(error &err, const int errNo=1)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
#define FatalErrorInFunction
const objectRegistry & db() const
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
void write(vtk::formatter &fmt, const Type &val, const label n=1)
virtual const labelUList & faceCells() const
Foam::fvPatchFieldMapper.
const fvPatch & patch() const
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Generic GeometricField class.
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Smooth ATC in cells next to a set of patches supplied by type.
virtual void autoMap(const fvPatchFieldMapper &)
const Time & time() const noexcept
const Boundary & boundaryField() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const volVectorField & U() const