This boundary condition provides a wall constraint on the kinematic turbulent thermal conductivity (i.e. alphat
) for atmospheric boundary layer modelling. It assumes a logarithmic distribution of the potential temperature within the first cell.
More...
Public Member Functions | |
TypeName ("atmAlphatkWallFunction") | |
Runtime type information. More... | |
atmAlphatkWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
atmAlphatkWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
atmAlphatkWallFunctionFvPatchScalarField (const atmAlphatkWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
atmAlphatkWallFunctionFvPatchScalarField (const atmAlphatkWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
virtual tmp< fvPatchScalarField > | clone () const |
Construct and return a clone. More... | |
atmAlphatkWallFunctionFvPatchScalarField (const atmAlphatkWallFunctionFvPatchScalarField &, 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 | updateCoeffs () |
Update the coefficients associated with the patch field. 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... | |
Protected Member Functions | |
virtual void | checkType () |
Check the type of the patch. More... | |
void | writeLocalEntries (Ostream &) const |
Write local wall function variables. More... | |
Protected Attributes | |
const scalar | Cmu_ |
Empirical model constant. More... | |
const scalar | kappa_ |
von Kármán constant More... | |
autoPtr< Function1< scalar > > | Pr_ |
Molecular Prandtl number. More... | |
autoPtr< PatchFunction1< scalar > > | Prt_ |
Turbulent Prandtl number field. More... | |
autoPtr< PatchFunction1< scalar > > | z0_ |
Surface roughness length [m]. More... | |
Static Protected Attributes | |
static scalar | tolerance_ = 0.01 |
Solution parameters. More... | |
static label | maxIters_ = 10 |
This boundary condition provides a wall constraint on the kinematic turbulent thermal conductivity (i.e. alphat
) for atmospheric boundary layer modelling. It assumes a logarithmic distribution of the potential temperature within the first cell.
Required fields:
alphat | Kinematic turbulent thermal conductivity [m2/s]
<patchName> { // Mandatory entries type atmAlphatkWallFunction; Pr <Function1<scalar>>; Prt <PatchFunction1<scalar>>; z0 <PatchFunction1<scalar>>; // Optional entries Cmu <scalar>; kappa <scalar>; // Inherited entries ... }
where the entries mean:
Property | Description | Type | Reqd | Deflt |
---|---|---|---|---|
type | Type name: atmAlphatkWallFunction | word | yes | - |
Pr | Molecular Prandtl number | Function1<scalar> | yes | - |
Prt | Turbulent Prandtl number | PatchFunction1<scalar> | yes | - |
z0 | Surface roughness length [m] | PatchFunction1<scalar> | yes | - |
Cmu | Empirical model constant | scalar | no | 0.09 |
kappa | von Kármán constant | scalar | no | 0.41 |
The inherited entries are elaborated in:
Definition at line 144 of file atmAlphatkWallFunctionFvPatchScalarField.H.
atmAlphatkWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 87 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References atmAlphatkWallFunctionFvPatchScalarField::checkType().
atmAlphatkWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Construct from patch, internal field and dictionary.
Definition at line 125 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References atmAlphatkWallFunctionFvPatchScalarField::checkType().
atmAlphatkWallFunctionFvPatchScalarField | ( | const atmAlphatkWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Construct by mapping given atmAlphatkWallFunctionFvPatchScalarField onto a new patch
Definition at line 105 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References atmAlphatkWallFunctionFvPatchScalarField::checkType().
Construct as copy.
Definition at line 160 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References atmAlphatkWallFunctionFvPatchScalarField::checkType().
atmAlphatkWallFunctionFvPatchScalarField | ( | const atmAlphatkWallFunctionFvPatchScalarField & | wfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct as copy setting internal field reference.
Definition at line 177 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References atmAlphatkWallFunctionFvPatchScalarField::checkType().
|
protectedvirtual |
Check the type of the patch.
Definition at line 48 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References Foam::abort(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, and Foam::nl.
Referenced by Foam::atan2(), atmAlphatkWallFunctionFvPatchScalarField::atmAlphatkWallFunctionFvPatchScalarField(), Foam::hypot(), Foam::max(), Foam::min(), Foam::operator+(), and Foam::operator-().
|
protected |
Write local wall function variables.
Definition at line 62 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References atmAlphatkWallFunctionFvPatchScalarField::Cmu_, atmAlphatkWallFunctionFvPatchScalarField::kappa_, os(), atmAlphatkWallFunctionFvPatchScalarField::Pr_, atmAlphatkWallFunctionFvPatchScalarField::Prt_, Ostream::writeEntryIfDifferent(), and atmAlphatkWallFunctionFvPatchScalarField::z0_.
Referenced by atmAlphatkWallFunctionFvPatchScalarField::write().
TypeName | ( | "atmAlphatkWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Construct and return a clone.
Definition at line 224 of file atmAlphatkWallFunctionFvPatchScalarField.H.
|
inlinevirtual |
Construct and return a clone setting internal field reference.
Definition at line 240 of file atmAlphatkWallFunctionFvPatchScalarField.H.
|
virtual |
Update the coefficients associated with the patch field.
Definition at line 197 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References atmAlphatkWallFunctionFvPatchScalarField::Cmu_, e, Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, IOobject::groupName(), k, atmAlphatkWallFunctionFvPatchScalarField::kappa_, Foam::log(), Foam::max(), Foam::pow025(), Pr(), atmAlphatkWallFunctionFvPatchScalarField::Pr_, turbulenceModel::propertiesName, Prt(), atmAlphatkWallFunctionFvPatchScalarField::Prt_, Foam::sqrt(), atmBoundaryLayerInletEpsilonFvPatchScalarField::updateCoeffs(), y, and atmAlphatkWallFunctionFvPatchScalarField::z0_.
|
virtual |
Map (and resize as needed) from self given a mapping object.
Definition at line 280 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References atmAlphatkWallFunctionFvPatchScalarField::Prt_, and atmAlphatkWallFunctionFvPatchScalarField::z0_.
|
virtual |
Reverse map the given fvPatchField onto this fvPatchField.
Definition at line 298 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References atmAlphatkWallFunctionFvPatchScalarField::Prt_, and atmAlphatkWallFunctionFvPatchScalarField::z0_.
|
virtual |
Write.
Definition at line 320 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References os(), ObukhovLength::write(), and atmAlphatkWallFunctionFvPatchScalarField::writeLocalEntries().
|
protected |
Empirical model constant.
Definition at line 153 of file atmAlphatkWallFunctionFvPatchScalarField.H.
Referenced by atmAlphatkWallFunctionFvPatchScalarField::updateCoeffs(), and atmAlphatkWallFunctionFvPatchScalarField::writeLocalEntries().
|
protected |
von Kármán constant
Definition at line 156 of file atmAlphatkWallFunctionFvPatchScalarField.H.
Referenced by atmAlphatkWallFunctionFvPatchScalarField::updateCoeffs(), and atmAlphatkWallFunctionFvPatchScalarField::writeLocalEntries().
Molecular Prandtl number.
Definition at line 159 of file atmAlphatkWallFunctionFvPatchScalarField.H.
Referenced by atmAlphatkWallFunctionFvPatchScalarField::updateCoeffs(), and atmAlphatkWallFunctionFvPatchScalarField::writeLocalEntries().
|
protected |
Turbulent Prandtl number field.
Definition at line 162 of file atmAlphatkWallFunctionFvPatchScalarField.H.
Referenced by atmAlphatkWallFunctionFvPatchScalarField::autoMap(), atmAlphatkWallFunctionFvPatchScalarField::rmap(), atmAlphatkWallFunctionFvPatchScalarField::updateCoeffs(), and atmAlphatkWallFunctionFvPatchScalarField::writeLocalEntries().
|
protected |
Surface roughness length [m].
Definition at line 165 of file atmAlphatkWallFunctionFvPatchScalarField.H.
Referenced by atmAlphatkWallFunctionFvPatchScalarField::autoMap(), atmAlphatkWallFunctionFvPatchScalarField::rmap(), atmAlphatkWallFunctionFvPatchScalarField::updateCoeffs(), and atmAlphatkWallFunctionFvPatchScalarField::writeLocalEntries().
|
staticprotected |
Solution parameters.
Definition at line 170 of file atmAlphatkWallFunctionFvPatchScalarField.H.
|
staticprotected |
Definition at line 171 of file atmAlphatkWallFunctionFvPatchScalarField.H.