This boundary condition provides a wall constraint on the turbulent kinematic viscosity, i.e. nut
, when using wall functions for rough walls, based on velocity, U
, using Spalding's law to give a continuous nut
profile to the wall (y+ = 0)
More...
Public Member Functions | |
TypeName ("nutUSpaldingWallFunction") | |
Runtime type information. More... | |
nutUSpaldingWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
nutUSpaldingWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
nutUSpaldingWallFunctionFvPatchScalarField (const nutUSpaldingWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
nutUSpaldingWallFunctionFvPatchScalarField (const nutUSpaldingWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
virtual tmp< fvPatchScalarField > | clone () const |
Construct and return a clone. More... | |
nutUSpaldingWallFunctionFvPatchScalarField (const nutUSpaldingWallFunctionFvPatchScalarField &, 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 | ~nutUSpaldingWallFunctionFvPatchScalarField () |
Destructor. More... | |
virtual tmp< scalarField > | yPlus () const |
Calculate and return the yPlus at the boundary. More... | |
virtual void | write (Ostream &os) const |
Write. More... | |
Public Member Functions inherited from nutWallFunctionFvPatchScalarField | |
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 y+ at the edge of the viscous sublayer. More... | |
virtual void | updateCoeffs () |
Update the coefficients associated with the patch field. More... | |
Protected Member Functions | |
virtual tmp< scalarField > | calcNut () const |
Uncomment in case of intrumentation. More... | |
virtual tmp< scalarField > | calcUTau (const scalarField &magGradU) const |
Calculate the friction velocity. More... | |
virtual tmp< scalarField > | calcUTau (const scalarField &magGradU, const label maxIter, scalarField &err) const |
virtual void | writeLocalEntries (Ostream &) const |
Write local wall function variables. More... | |
Protected Member Functions inherited from nutWallFunctionFvPatchScalarField | |
virtual const volVectorField & | U (const turbulenceModel &turb) const |
virtual void | checkType () |
Check the type of the patch. More... | |
Protected Attributes | |
const label | maxIter_ |
Max iterations in calcNut. More... | |
const scalar | tolerance_ |
Convergence tolerance. More... | |
Protected Attributes inherited from nutWallFunctionFvPatchScalarField | |
word | UName_ |
Name of velocity field. More... | |
scalar | Cmu_ |
Cmu coefficient. More... | |
scalar | kappa_ |
Von Karman constant. More... | |
scalar | E_ |
E coefficient. More... | |
scalar | yPlusLam_ |
Estimated y+ value at the edge of the viscous sublayer. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from nutWallFunctionFvPatchScalarField | |
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) |
Calculate the y+ at the edge of the viscous sublayer. More... | |
This boundary condition provides a wall constraint on the turbulent kinematic viscosity, i.e. nut
, when using wall functions for rough walls, based on velocity, U
, using Spalding's law to give a continuous nut
profile to the wall (y+ = 0)
\[ y^+ = u^+ + \frac{1}{E} \left[exp(\kappa u^+) - 1 - \kappa u^+\, - 0.5 (\kappa u^+)^2 - \frac{1}{6} (\kappa u^+)^3\right] \]
where
\( y^+ \) | = | non-dimensional position |
\( u^+ \) | = | non-dimensional velocity |
\( \kappa \) | = | Von Karman constant |
Property | Description | Required | Default value |
---|---|---|---|
Cmu | Model coefficient | no | 0.09 |
kappa | Von Karman constant | no | 0.41 |
E | Model coefficient | no | 9.8 |
maxIter | Number of Newton-Raphson iterations | no | 10 |
tolerance | Convergence tolerance | no | 0.0001 |
Example of the boundary condition specification:
<patchName> { // Mandatory entries type nutUSpaldingWallFunction; // Optional entries }
This can be avoided by overriding the tolerance. This also switches on a pre-detection whether the current nut already satisfies the turbulence conditions and if so does not change it at all. This means that the nut only changes if it 'has really changed'. This probably should be used with a tight tolerance, e.g.
maxIter 100; tolerance 1e-7;
to make sure to kick every iteration.
Definition at line 152 of file nutUSpaldingWallFunctionFvPatchScalarField.H.
nutUSpaldingWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 209 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
Referenced by nutUSpaldingWallFunctionFvPatchScalarField::clone().
nutUSpaldingWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Construct from patch, internal field and dictionary.
Definition at line 245 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
nutUSpaldingWallFunctionFvPatchScalarField | ( | const nutUSpaldingWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Construct by mapping given nutUSpaldingWallFunctionFvPatchScalarField onto a new patch
Definition at line 226 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
nutUSpaldingWallFunctionFvPatchScalarField | ( | const nutUSpaldingWallFunctionFvPatchScalarField & | wfpsf | ) |
Construct as copy.
Definition at line 263 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
nutUSpaldingWallFunctionFvPatchScalarField | ( | const nutUSpaldingWallFunctionFvPatchScalarField & | wfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct as copy setting internal field reference.
Definition at line 279 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
|
virtual |
Destructor.
Definition at line 297 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
|
protectedvirtual |
Uncomment in case of intrumentation.
mutable uint64_t invocations_; mutable uint64_t nontrivial_; mutable uint64_t nonconvergence_; mutable uint64_t iterations_; Calculate the turbulence viscosity
Implements nutWallFunctionFvPatchScalarField.
Definition at line 39 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), nutUSpaldingWallFunctionFvPatchScalarField::calcUTau(), forAll, Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), Foam::max(), turbulenceModel::nu(), nutWallFunctionFvPatchScalarField::nutw(), Foam::foamVersion::patch, turbulenceModel::propertiesName, tmp< T >::ref(), Foam::sqr(), nutUSpaldingWallFunctionFvPatchScalarField::tolerance_, and nutWallFunctionFvPatchScalarField::U().
|
protectedvirtual |
Calculate the friction velocity.
Definition at line 94 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
Referenced by nutUSpaldingWallFunctionFvPatchScalarField::calcNut().
|
protectedvirtual |
Calculate the friction velocity and number of iterations for convergence
Definition at line 105 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
References Foam::exp(), f(), forAll, Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), magUp, Foam::max(), Foam::min(), turbulenceModel::nu(), Foam::foamVersion::patch, fvPatchField< Type >::patchInternalField(), turbulenceModel::propertiesName, tmp< T >::ref(), Foam::sqr(), Foam::sqrt(), U, uTau, y, turbulenceModel::y(), and Foam::Zero.
|
protectedvirtual |
Write local wall function variables.
Reimplemented from nutWallFunctionFvPatchScalarField.
Definition at line 194 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
References Ostream::writeEntryIfDifferent(), and nutWallFunctionFvPatchScalarField::writeLocalEntries().
TypeName | ( | "nutUSpaldingWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Construct and return a clone.
Definition at line 235 of file nutUSpaldingWallFunctionFvPatchScalarField.H.
References nutUSpaldingWallFunctionFvPatchScalarField::nutUSpaldingWallFunctionFvPatchScalarField().
|
inlinevirtual |
Construct and return a clone setting internal field reference.
Definition at line 252 of file nutUSpaldingWallFunctionFvPatchScalarField.H.
References nutUSpaldingWallFunctionFvPatchScalarField::nutUSpaldingWallFunctionFvPatchScalarField().
|
virtual |
Calculate and return the yPlus at the boundary.
Implements nutWallFunctionFvPatchScalarField.
Definition at line 318 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
References Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), turbulenceModel::nu(), Foam::foamVersion::patch, turbulenceModel::propertiesName, fvPatchField< Type >::snGrad(), U, y, and turbulenceModel::y().
|
virtual |
Write.
Reimplemented from nutWallFunctionFvPatchScalarField.
Definition at line 340 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
References fvPatchField< Type >::write().
|
protected |
Max iterations in calcNut.
Definition at line 161 of file nutUSpaldingWallFunctionFvPatchScalarField.H.
|
protected |
Convergence tolerance.
Definition at line 164 of file nutUSpaldingWallFunctionFvPatchScalarField.H.
Referenced by nutUSpaldingWallFunctionFvPatchScalarField::calcNut().