This boundary condition provides a wall constraint on the turbulent viscosity (i.e. nut
) based on the turbulent kinetic energy (i.e. k
) for atmospheric boundary layer modelling. It is designed to be used in conjunction with the atmBoundaryLayerInletVelocity
boundary condition.
More...
Public Member Functions | |
TypeName ("atmNutkWallFunction") | |
Runtime type information. More... | |
atmNutkWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
atmNutkWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
atmNutkWallFunctionFvPatchScalarField (const atmNutkWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
atmNutkWallFunctionFvPatchScalarField (const atmNutkWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
virtual tmp< fvPatchScalarField > | clone () const |
Construct and return a clone. More... | |
atmNutkWallFunctionFvPatchScalarField (const atmNutkWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
Construct as copy setting internal field reference. More... | |
virtual tmp< fvPatchScalarField > | clone (const DimensionedField< scalar, volMesh > &iF) const |
Construct and return a clone setting internal field reference. More... | |
virtual void | autoMap (const fvPatchFieldMapper &) |
Map (and resize as needed) from self given a mapping object. More... | |
virtual void | rmap (const fvPatchScalarField &, const labelList &) |
Reverse map the given fvPatchField onto this fvPatchField. More... | |
virtual void | write (Ostream &) const |
Write. More... | |
![]() | |
TypeName ("nutkWallFunction") | |
Runtime type information. More... | |
nutkWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
nutkWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
nutkWallFunctionFvPatchScalarField (const nutkWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
nutkWallFunctionFvPatchScalarField (const nutkWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
nutkWallFunctionFvPatchScalarField (const nutkWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
Construct as copy setting internal field reference. More... | |
virtual tmp< scalarField > | yPlus () const |
Calculate and return the yPlus at the boundary. More... | |
![]() | |
TypeName ("nutWallFunction") | |
Runtime type information. More... | |
nutWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
nutWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
nutWallFunctionFvPatchScalarField (const nutWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
Construct as copy setting internal field reference. More... | |
scalar | Cmu () const |
Return Cmu. More... | |
scalar | kappa () const |
Return kappa. More... | |
scalar | E () const |
Return E. More... | |
scalar | yPlusLam () const |
Return the estimated y+ at the two-sublayer intersection. More... | |
scalar | blend (const scalar nutVis, const scalar nutLog, const scalar yPlus) const |
Return the blended nut according to the chosen blending treatment. More... | |
virtual void | updateCoeffs () |
Update the coefficients associated with the patch field. More... | |
Protected Member Functions | |
virtual tmp< scalarField > | calcNut () const |
Calculate the turbulent viscosity. More... | |
![]() | |
virtual const volVectorField & | U (const turbulenceModel &turb) const |
virtual void | checkType () |
Check the type of the patch. More... | |
virtual void | writeLocalEntries (Ostream &) const |
Write local wall function variables. More... | |
Additional Inherited Members | |
![]() | |
static const nutWallFunctionFvPatchScalarField & | nutw (const turbulenceModel &turbModel, const label patchi) |
Return the nut patchField for the given wall patch. More... | |
static scalar | yPlusLam (const scalar kappa, const scalar E) |
Estimate the y+ at the intersection of the two sublayers. More... | |
![]() | |
enum | blendingType { STEPWISE, MAX, BINOMIAL, EXPONENTIAL } |
Options for the blending treatment of viscous and inertial sublayers. More... | |
![]() | |
enum blendingType | blending_ |
Blending treatment (default = blendingType::STEPWISE) More... | |
const scalar | n_ |
word | UName_ |
Name of velocity field. More... | |
scalar | Cmu_ |
Empirical model coefficient. More... | |
scalar | kappa_ |
von Kármán constant More... | |
scalar | E_ |
Wall roughness parameter. More... | |
scalar | yPlusLam_ |
![]() | |
static const Enum< blendingType > | blendingTypeNames |
Names for blendingType. More... | |
This boundary condition provides a wall constraint on the turbulent viscosity (i.e. nut
) based on the turbulent kinetic energy (i.e. k
) for atmospheric boundary layer modelling. It is designed to be used in conjunction with the atmBoundaryLayerInletVelocity
boundary condition.
The governing equation of the boundary condition:
\[ u = \frac{u^*}{\kappa} ln \left(\frac{z + z_0}{z_0}\right) \]
where
\( u^* \) | = | Friction velocity |
\( \kappa \) | = | von Kármán constant |
\( z_0 \) | = | Surface roughness length [m] |
\( z \) | = | Ground-normal coordinate |
Required fields:
nut | Turbulent viscosity [m2/s] k | Turbulent kinetic energy [m2/s2]
References:
Theoretical expressions (tag:HW): Hargreaves, D. M., & Wright, N. G. (2007). On the use of the k–ε model in commercial CFD software to model the neutral atmospheric boundary layer. J. of wind engineering and industrial aerodynamics, 95(5), 355-369. DOI:10.1016/j.jweia.2006.08.002
Required fields:
nut | Turbulent viscosity [m2/s] k | Turbulent kinetic energy [m2/s2]
<patchName> { // Mandatory entries (unmodifiable) type atmNutkWallFunction; // Mandatory entries (runtime modifiable) z0 uniform 0.001; // Optional entries (unmodifiable) boundNut false; // Optional (inherited) entries ... }
where the entries mean:
Property | Description | Type | Reqd | Dflt |
---|---|---|---|---|
type | Type name: atmNutkWallFunction | word | yes | - |
z0 | Surface roughness length [m] | PatchFunction1<scalar> | yes | - |
boundNut | Flag to zero-bound nut near wall | bool | no | false |
The inherited entries are elaborated in:
boundNut
entry is set false
for backward compatibility reasons.nutkAtmRoughWallFunction
was renamed to atmNutkWallFunction
.Definition at line 173 of file atmNutkWallFunctionFvPatchScalarField.H.
atmNutkWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 112 of file atmNutkWallFunctionFvPatchScalarField.C.
Referenced by atmNutkWallFunctionFvPatchScalarField::clone().
atmNutkWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Construct from patch, internal field and dictionary.
Definition at line 138 of file atmNutkWallFunctionFvPatchScalarField.C.
atmNutkWallFunctionFvPatchScalarField | ( | const atmNutkWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Construct by mapping given atmNutkWallFunctionFvPatchScalarField onto a new patch
Definition at line 124 of file atmNutkWallFunctionFvPatchScalarField.C.
atmNutkWallFunctionFvPatchScalarField | ( | const atmNutkWallFunctionFvPatchScalarField & | rwfpsf | ) |
Construct as copy.
Definition at line 151 of file atmNutkWallFunctionFvPatchScalarField.C.
atmNutkWallFunctionFvPatchScalarField | ( | const atmNutkWallFunctionFvPatchScalarField & | rwfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct as copy setting internal field reference.
Definition at line 162 of file atmNutkWallFunctionFvPatchScalarField.C.
|
protectedvirtual |
Calculate the turbulent viscosity.
Reimplemented from nutkWallFunctionFvPatchScalarField.
Definition at line 44 of file atmNutkWallFunctionFvPatchScalarField.C.
References nutWallFunctionFvPatchScalarField::Cmu_, Foam::constant::electromagnetic::e, Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::constant::atomic::group, IOobject::groupName(), k, nutWallFunctionFvPatchScalarField::kappa_, Foam::log(), Foam::max(), tmp< T >::New(), nutWallFunctionFvPatchScalarField::nutw(), Foam::foamVersion::patch, Foam::pow025(), turbulenceModel::propertiesName, Foam::sqrt(), y, and nutkWallFunctionFvPatchScalarField::yPlus().
TypeName | ( | "atmNutkWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Construct and return a clone.
Reimplemented from nutkWallFunctionFvPatchScalarField.
Definition at line 236 of file atmNutkWallFunctionFvPatchScalarField.H.
References atmNutkWallFunctionFvPatchScalarField::atmNutkWallFunctionFvPatchScalarField().
|
inlinevirtual |
Construct and return a clone setting internal field reference.
Reimplemented from nutkWallFunctionFvPatchScalarField.
Definition at line 253 of file atmNutkWallFunctionFvPatchScalarField.H.
References atmNutkWallFunctionFvPatchScalarField::atmNutkWallFunctionFvPatchScalarField().
|
virtual |
Map (and resize as needed) from self given a mapping object.
Definition at line 176 of file atmNutkWallFunctionFvPatchScalarField.C.
|
virtual |
Reverse map the given fvPatchField onto this fvPatchField.
Definition at line 186 of file atmNutkWallFunctionFvPatchScalarField.C.
|
virtual |
Write.
Reimplemented from nutWallFunctionFvPatchScalarField.
Definition at line 200 of file atmNutkWallFunctionFvPatchScalarField.C.
References os(), fvPatchField< Type >::write(), and nutWallFunctionFvPatchScalarField::writeLocalEntries().