This boundary condition provides a kinematic turbulent thermal conductivity for using wall functions, using the Jayatilleke 'P' function. More...
Public Member Functions | |
TypeName ("alphatJayatillekeWallFunction") | |
Runtime type information. More... | |
alphatJayatillekeWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
alphatJayatillekeWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
alphatJayatillekeWallFunctionFvPatchScalarField (const alphatJayatillekeWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
Construct by mapping given alphatJayatillekeWallFunctionFvPatchScalarField onto a new patch. More... | |
alphatJayatillekeWallFunctionFvPatchScalarField (const alphatJayatillekeWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
alphatJayatillekeWallFunctionFvPatchScalarField (const alphatJayatillekeWallFunctionFvPatchScalarField &, 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 | write (Ostream &) const |
Write. More... | |
Protected Member Functions | |
virtual void | checkType () |
Check the type of the patch. More... | |
tmp< scalarField > | yPlus (const turbulenceModel &turbModel) const |
Return the patch y+. More... | |
scalar | Psmooth (const scalar Prat) const |
`P' function More... | |
scalar | yPlusTherm (const scalar P, const scalar Prat) const |
Calculate y+ at the edge of the thermal laminar sublayer. More... | |
Protected Attributes | |
scalar | Prt_ |
Turbulent Prandtl number. More... | |
scalar | kappa_ |
Von Karman constant. More... | |
scalar | E_ |
E coefficient. More... | |
Static Protected Attributes | |
static scalar | tolerance_ = 0.01 |
static label | maxIters_ = 10 |
This boundary condition provides a kinematic turbulent thermal conductivity for using wall functions, using the Jayatilleke 'P' function.
<patchName> { // Mandatory entries type alphatJayatillekeWallFunction; // Optional entries Prt <scalar>; kappa <scalar>; E <scalar>; // Inherited entries ... }
where the entries mean:
Property | Description | Type | Reqd | Deflt |
---|---|---|---|---|
type | Type name: alphatJayatillekeWallFunction | word | yes | - |
Prt | Turbulent Prandtl number | scalar | no | 0.85 |
kappa | von Karman constant | scalar | no | 0.41 |
E | Wall roughness parameter | scalar | no | 9.8 |
The inherited entries are elaborated in:
Definition at line 117 of file alphatJayatillekeWallFunctionFvPatchScalarField.H.
alphatJayatillekeWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 133 of file alphatJayatillekeWallFunctionFvPatchScalarField.C.
alphatJayatillekeWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Construct from patch, internal field and dictionary.
Definition at line 167 of file alphatJayatillekeWallFunctionFvPatchScalarField.C.
alphatJayatillekeWallFunctionFvPatchScalarField | ( | const alphatJayatillekeWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Construct by mapping given alphatJayatillekeWallFunctionFvPatchScalarField onto a new patch.
Definition at line 149 of file alphatJayatillekeWallFunctionFvPatchScalarField.C.
alphatJayatillekeWallFunctionFvPatchScalarField | ( | const alphatJayatillekeWallFunctionFvPatchScalarField & | wfpsf | ) |
Construct as copy.
Definition at line 184 of file alphatJayatillekeWallFunctionFvPatchScalarField.C.
alphatJayatillekeWallFunctionFvPatchScalarField | ( | const alphatJayatillekeWallFunctionFvPatchScalarField & | wfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct as copy setting internal field reference.
Definition at line 199 of file alphatJayatillekeWallFunctionFvPatchScalarField.C.
|
protectedvirtual |
Check the type of the patch.
Definition at line 43 of file alphatJayatillekeWallFunctionFvPatchScalarField.C.
References Foam::abort(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, Foam::nl, and Foam::foamVersion::patch.
|
protected |
Return the patch y+.
Definition at line 58 of file alphatJayatillekeWallFunctionFvPatchScalarField.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), Foam::mag(), TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >::nu(), turbulenceModel::nuEff(), nut, turbulenceModel::nut(), Foam::foamVersion::patch, fvPatchField< Type >::snGrad(), Foam::sqrt(), turbulenceModel::U(), y, turbulenceModel::y(), and nutWallFunctionFvPatchScalarField::yPlus().
Referenced by alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs().
|
protected |
`P' function
Definition at line 89 of file alphatJayatillekeWallFunctionFvPatchScalarField.C.
References Foam::exp(), and Foam::pow().
Referenced by alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs().
|
protected |
Calculate y+ at the edge of the thermal laminar sublayer.
Definition at line 98 of file alphatJayatillekeWallFunctionFvPatchScalarField.C.
References f(), Foam::log(), and Foam::mag().
Referenced by alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs().
TypeName | ( | "alphatJayatillekeWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Return a clone.
Definition at line 236 of file alphatJayatillekeWallFunctionFvPatchScalarField.H.
References fvPatchField< Type >::Clone().
|
inlinevirtual |
Clone with an internal field reference.
Definition at line 245 of file alphatJayatillekeWallFunctionFvPatchScalarField.H.
References fvPatchField< Type >::Clone().
|
virtual |
Update the coefficients associated with the patch field.
Definition at line 215 of file alphatJayatillekeWallFunctionFvPatchScalarField.C.
References Foam::dimless, alphatJayatillekeWallFunctionFvPatchScalarField::E_, forAll, Foam::constant::atomic::group, IOobject::groupName(), alphatJayatillekeWallFunctionFvPatchScalarField::kappa_, Foam::log(), Foam::max(), TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >::nu(), nu, Foam::foamVersion::patch, Pr(), turbulenceModel::propertiesName, alphatJayatillekeWallFunctionFvPatchScalarField::Prt_, alphatJayatillekeWallFunctionFvPatchScalarField::Psmooth(), transportProperties(), fvPatchField< Type >::updateCoeffs(), alphatJayatillekeWallFunctionFvPatchScalarField::yPlus(), and alphatJayatillekeWallFunctionFvPatchScalarField::yPlusTherm().
|
virtual |
Write.
Definition at line 282 of file alphatJayatillekeWallFunctionFvPatchScalarField.C.
References alphatJayatillekeWallFunctionFvPatchScalarField::E_, alphatJayatillekeWallFunctionFvPatchScalarField::kappa_, os(), alphatJayatillekeWallFunctionFvPatchScalarField::Prt_, fvPatchField< Type >::write(), and fvPatchField< Type >::writeValueEntry().
|
protected |
Turbulent Prandtl number.
Definition at line 128 of file alphatJayatillekeWallFunctionFvPatchScalarField.H.
Referenced by alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs(), and alphatJayatillekeWallFunctionFvPatchScalarField::write().
|
protected |
Von Karman constant.
Definition at line 133 of file alphatJayatillekeWallFunctionFvPatchScalarField.H.
Referenced by alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs(), and alphatJayatillekeWallFunctionFvPatchScalarField::write().
|
protected |
E coefficient.
Definition at line 138 of file alphatJayatillekeWallFunctionFvPatchScalarField.H.
Referenced by alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs(), and alphatJayatillekeWallFunctionFvPatchScalarField::write().
|
staticprotected |
Definition at line 143 of file alphatJayatillekeWallFunctionFvPatchScalarField.H.
|
staticprotected |
Definition at line 144 of file alphatJayatillekeWallFunctionFvPatchScalarField.H.