This boundary condition provides a wall constraint on the turbulent viscosity (i.e. nut
) based on the turbulent kinetic energy (i.e. k
) and velocity (i.e. U
) for atmospheric boundary layer modelling.
More...
Public Member Functions | |
TypeName ("atmNutWallFunction") | |
Runtime type information. More... | |
atmNutWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
atmNutWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
atmNutWallFunctionFvPatchScalarField (const atmNutWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
atmNutWallFunctionFvPatchScalarField (const atmNutWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
virtual tmp< fvPatchScalarField > | clone () const |
Construct and return a clone. More... | |
atmNutWallFunctionFvPatchScalarField (const atmNutWallFunctionFvPatchScalarField &, 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
) and velocity (i.e. U
) for atmospheric boundary layer modelling.
The governing equation of the boundary condition:
\[ \tau_w = {U^*_u} {U^*_k} \]
with
\[ {U^*_u} = \frac{\kappa U_w}{ln(z_p / z_0)} \]
\[ {U^*_k} = C_{\mu}^{1/4} \sqrt{k} \]
where
\( \tau_w \) | = | wall shear stress |
\( U^*_u \) | = | local friction velocity based on near-ground velocity |
\( U^*_k \) | = | local friction velocity based on near-ground k |
\( \kappa \) | = | von Kármán constant |
\( U_w \) | = | near-ground velocity |
\( z_p \) | = | vertical coordinate |
\( z_0 \) | = | surface roughness length [m] |
\( C_mu \) | = | empirical model constant |
\( k \) | = | turbulent kinetic energy |
References:
Theoretical expressions (tags:RH, SBJM, SM): Richards, P. J., & Hoxey, R. P. (1993). Appropriate boundary conditions for computational wind engineering models using the k-ε turbulence model. In Computational Wind Engineering 1 (pp. 145-153). DOI:10.1016/B978-0-444-81688-7.50018-8 Sørensen, N. N., Bechmann, A., Johansen, J., Myllerup, L., Botha, P., Vinther, S., & Nielsen, B. S. (2007). Identification of severe wind conditions using a Reynolds Averaged Navier-Stokes solver. In Journal of Physics: Conference series (Vol. 75, No. 1, p. 012053). DOI:10.1088/1742-6596/75/1/012053 Sumner, J., & Masson, C. (2012). k−ε simulations of the neutral atmospheric boundary layer: analysis and correction of discretization errors on practical grids. International journal for numerical methods in fluids, 70(6), 724-741. DOI:10.1002/fld.2709
Required fields:
nut | Turbulent viscosity [m2/s] k | Turbulent kinetic energy [m2/s2]
<patchName> { // Mandatory entries (unmodifiable) type atmNutWallFunction; z0Min 0.001; // Mandatory entries (runtime modifiable) z0 uniform 0.001; // Optional (inherited) entries ... }
where the entries mean:
Property | Description | Type | Reqd | Dflt |
---|---|---|---|---|
type | Type name: nutAtmWallFunction | word | yes | - |
z0Min | Minimum surface roughness length [m] | scalar | yes | - |
z0 | Surface roughness length [m] | PatchFunction1<scalar> | yes | - |
The inherited entries are elaborated in:
Definition at line 211 of file atmNutWallFunctionFvPatchScalarField.H.
atmNutWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 117 of file atmNutWallFunctionFvPatchScalarField.C.
Referenced by atmNutWallFunctionFvPatchScalarField::clone().
atmNutWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Construct from patch, internal field and dictionary.
Definition at line 143 of file atmNutWallFunctionFvPatchScalarField.C.
atmNutWallFunctionFvPatchScalarField | ( | const atmNutWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Construct by mapping given atmNutWallFunctionFvPatchScalarField onto a new patch
Definition at line 129 of file atmNutWallFunctionFvPatchScalarField.C.
atmNutWallFunctionFvPatchScalarField | ( | const atmNutWallFunctionFvPatchScalarField & | rwfpsf | ) |
Construct as copy.
Definition at line 164 of file atmNutWallFunctionFvPatchScalarField.C.
atmNutWallFunctionFvPatchScalarField | ( | const atmNutWallFunctionFvPatchScalarField & | rwfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct as copy setting internal field reference.
Definition at line 175 of file atmNutWallFunctionFvPatchScalarField.C.
|
protectedvirtual |
Calculate the turbulent viscosity.
Reimplemented from nutkWallFunctionFvPatchScalarField.
Definition at line 43 of file atmNutWallFunctionFvPatchScalarField.C.
References nutWallFunctionFvPatchScalarField::Cmu_, Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::constant::atomic::group, IOobject::groupName(), k, nutWallFunctionFvPatchScalarField::kappa_, Foam::log(), Foam::mag(), magUp, Foam::max(), tmp< T >::New(), nutWallFunctionFvPatchScalarField::nutw(), Foam::foamVersion::patch, fvPatchField< Type >::patchInternalField(), Foam::pow025(), turbulenceModel::propertiesName, Foam::sqrt(), tauw, and y.
TypeName | ( | "atmNutWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Construct and return a clone.
Reimplemented from nutkWallFunctionFvPatchScalarField.
Definition at line 274 of file atmNutWallFunctionFvPatchScalarField.H.
References atmNutWallFunctionFvPatchScalarField::atmNutWallFunctionFvPatchScalarField().
|
inlinevirtual |
Construct and return a clone setting internal field reference.
Reimplemented from nutkWallFunctionFvPatchScalarField.
Definition at line 291 of file atmNutWallFunctionFvPatchScalarField.H.
References atmNutWallFunctionFvPatchScalarField::atmNutWallFunctionFvPatchScalarField().
|
virtual |
Map (and resize as needed) from self given a mapping object.
Definition at line 189 of file atmNutWallFunctionFvPatchScalarField.C.
|
virtual |
Reverse map the given fvPatchField onto this fvPatchField.
Definition at line 199 of file atmNutWallFunctionFvPatchScalarField.C.
|
virtual |
Write.
Reimplemented from nutWallFunctionFvPatchScalarField.
Definition at line 213 of file atmNutWallFunctionFvPatchScalarField.C.
References os(), fvPatchField< Type >::write(), and nutWallFunctionFvPatchScalarField::writeLocalEntries().