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 &) | |
Construct by mapping given atmAlphatkWallFunctionFvPatchScalarField onto a new patch. More... | |
atmAlphatkWallFunctionFvPatchScalarField (const atmAlphatkWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
atmAlphatkWallFunctionFvPatchScalarField (const atmAlphatkWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
Construct as copy setting internal field reference. More... | |
virtual tmp< fvPatchField< scalar > > | clone () const |
Return a clone. More... | |
virtual tmp< fvPatchField< scalar > > | clone (const DimensionedField< scalar, volMesh > &iF) const |
Clone with an 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.
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 140 of file atmAlphatkWallFunctionFvPatchScalarField.H.
atmAlphatkWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 82 of file atmAlphatkWallFunctionFvPatchScalarField.C.
atmAlphatkWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Construct from patch, internal field and dictionary.
Definition at line 120 of file atmAlphatkWallFunctionFvPatchScalarField.C.
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 100 of file atmAlphatkWallFunctionFvPatchScalarField.C.
Construct as copy.
Definition at line 155 of file atmAlphatkWallFunctionFvPatchScalarField.C.
atmAlphatkWallFunctionFvPatchScalarField | ( | const atmAlphatkWallFunctionFvPatchScalarField & | wfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct as copy setting internal field reference.
Definition at line 172 of file atmAlphatkWallFunctionFvPatchScalarField.C.
|
protectedvirtual |
Check the type of the patch.
Definition at line 41 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References Foam::abort(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, Foam::nl, and Foam::foamVersion::patch.
|
protected |
Write local wall function variables.
Definition at line 56 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References os(), and Ostream::writeEntryIfDifferent().
Referenced by atmAlphatkWallFunctionFvPatchScalarField::write().
TypeName | ( | "atmAlphatkWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Return a clone.
Definition at line 257 of file atmAlphatkWallFunctionFvPatchScalarField.H.
References fvPatchField< Type >::Clone().
|
inlinevirtual |
Clone with an internal field reference.
Definition at line 266 of file atmAlphatkWallFunctionFvPatchScalarField.H.
References fvPatchField< Type >::Clone().
|
virtual |
Update the coefficients associated with the patch field.
Definition at line 190 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References atmAlphatkWallFunctionFvPatchScalarField::Cmu_, Foam::constant::electromagnetic::e, Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::constant::atomic::group, IOobject::groupName(), k, atmAlphatkWallFunctionFvPatchScalarField::kappa_, Foam::log(), Foam::max(), Foam::foamVersion::patch, Foam::pow025(), Pr(), atmAlphatkWallFunctionFvPatchScalarField::Pr_, turbulenceModel::propertiesName, Prt(), atmAlphatkWallFunctionFvPatchScalarField::Prt_, Foam::sqrt(), fvPatchField< Type >::updateCoeffs(), y, and atmAlphatkWallFunctionFvPatchScalarField::z0_.
|
virtual |
Map (and resize as needed) from self given a mapping object.
Definition at line 274 of file atmAlphatkWallFunctionFvPatchScalarField.C.
|
virtual |
Reverse map the given fvPatchField onto this fvPatchField.
Definition at line 292 of file atmAlphatkWallFunctionFvPatchScalarField.C.
|
virtual |
Write.
Definition at line 313 of file atmAlphatkWallFunctionFvPatchScalarField.C.
References os(), fvPatchField< Type >::write(), atmAlphatkWallFunctionFvPatchScalarField::writeLocalEntries(), and fvPatchField< Type >::writeValueEntry().
|
protected |
Empirical model constant.
Definition at line 151 of file atmAlphatkWallFunctionFvPatchScalarField.H.
Referenced by atmAlphatkWallFunctionFvPatchScalarField::updateCoeffs().
|
protected |
von Kármán constant
Definition at line 156 of file atmAlphatkWallFunctionFvPatchScalarField.H.
Referenced by atmAlphatkWallFunctionFvPatchScalarField::updateCoeffs().
Molecular Prandtl number.
Definition at line 161 of file atmAlphatkWallFunctionFvPatchScalarField.H.
Referenced by atmAlphatkWallFunctionFvPatchScalarField::updateCoeffs().
|
protected |
Turbulent Prandtl number field.
Definition at line 166 of file atmAlphatkWallFunctionFvPatchScalarField.H.
Referenced by atmAlphatkWallFunctionFvPatchScalarField::updateCoeffs().
|
protected |
Surface roughness length [m].
Definition at line 171 of file atmAlphatkWallFunctionFvPatchScalarField.H.
Referenced by atmAlphatkWallFunctionFvPatchScalarField::updateCoeffs().
|
staticprotected |
Solution parameters.
Definition at line 178 of file atmAlphatkWallFunctionFvPatchScalarField.H.
|
staticprotected |
Definition at line 179 of file atmAlphatkWallFunctionFvPatchScalarField.H.