This boundary condition provides a wall constraint on the turbulent viscosity, i.e. nut
, based on velocity, i.e. U
. Using Spalding's law gives a continuous nut
profile to the wall.
More...
Protected Member Functions | |
virtual tmp< scalarField > | calcNut () const |
virtual tmp< scalarField > | calcUTau (const scalarField &magGradU) const |
virtual tmp< scalarField > | calcUTau (const scalarField &magGradU, const label maxIter, scalarField &err) const |
virtual void | writeLocalEntries (Ostream &) const |
![]() | |
virtual const volVectorField & | U (const turbulenceModel &turb) const |
virtual void | checkType () |
Protected Attributes | |
const label | maxIter_ |
const scalar | tolerance_ |
![]() | |
enum blendingType | blending_ |
const scalar | n_ |
word | UName_ |
scalar | Cmu_ |
scalar | kappa_ |
scalar | E_ |
scalar | yPlusLam_ |
Additional Inherited Members | |
![]() | |
static const nutWallFunctionFvPatchScalarField & | nutw (const turbulenceModel &turbModel, const label patchi) |
static scalar | yPlusLam (const scalar kappa, const scalar E) |
![]() | |
enum | blendingType { STEPWISE, MAX, BINOMIAL, EXPONENTIAL } |
![]() | |
static const Enum< blendingType > | blendingTypeNames |
This boundary condition provides a wall constraint on the turbulent viscosity, i.e. nut
, based on velocity, i.e. U
. Using Spalding's law gives a continuous nut
profile to the wall.
where
![]() | = | Non-dimensional position |
![]() | = | Non-dimensional velocity |
![]() | = | von Kármán constant |
<patchName> { // Mandatory entries (unmodifiable) type nutUSpaldingWallFunction; // Optional entries (unmodifiable) maxIter 10; tolerance 0.0001; // Optional (inherited) entries ... }
where the entries mean:
Property | Description | Type | Req'd | Dflt |
---|---|---|---|---|
type | Type name: nutUBlendedWallFunction | word | yes | - |
maxIter | Number of Newton-Raphson iterations | label | no | 10 |
tolerance | Convergence tolerance | scalar | no | 0.0001 |
The inherited entries are elaborated in:
correctNut()
(called through turbulence->validate
) returns a slightly different value every time it is called. This is since the seed for the Newton-Raphson iteration uses the current value of *this
(=nut
).Definition at line 141 of file nutUSpaldingWallFunctionFvPatchScalarField.H.
nutUSpaldingWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Definition at line 202 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
Referenced by nutUSpaldingWallFunctionFvPatchScalarField::clone().
nutUSpaldingWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Definition at line 238 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
nutUSpaldingWallFunctionFvPatchScalarField | ( | const nutUSpaldingWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Definition at line 219 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
nutUSpaldingWallFunctionFvPatchScalarField | ( | const nutUSpaldingWallFunctionFvPatchScalarField & | wfpsf | ) |
Definition at line 256 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
nutUSpaldingWallFunctionFvPatchScalarField | ( | const nutUSpaldingWallFunctionFvPatchScalarField & | wfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Definition at line 272 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
|
virtual |
Definition at line 290 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
|
protectedvirtual |
Implements nutWallFunctionFvPatchScalarField.
Definition at line 32 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
References GeometricField::boundaryField(), nutUSpaldingWallFunctionFvPatchScalarField::calcUTau(), forAll, Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), Foam::max(), turbulenceModel::nu(), nutWallFunctionFvPatchScalarField::nutw(), Foam::foamVersion::patch, turbulenceModel::propertiesName, tmp::ref(), Foam::sqr(), nutUSpaldingWallFunctionFvPatchScalarField::tolerance_, and nutWallFunctionFvPatchScalarField::U().
|
protectedvirtual |
Definition at line 87 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
Referenced by nutUSpaldingWallFunctionFvPatchScalarField::calcNut().
|
protectedvirtual |
Definition at line 98 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::patchInternalField(), turbulenceModel::propertiesName, tmp::ref(), Foam::sqr(), Foam::sqrt(), U, uTau, y, turbulenceModel::y(), and Foam::Zero.
|
protectedvirtual |
Reimplemented from nutWallFunctionFvPatchScalarField.
Definition at line 187 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
References os(), Ostream::writeEntryIfDifferent(), and nutWallFunctionFvPatchScalarField::writeLocalEntries().
TypeName | ( | "nutUSpaldingWallFunction" | ) |
|
inlinevirtual |
Definition at line 224 of file nutUSpaldingWallFunctionFvPatchScalarField.H.
References nutUSpaldingWallFunctionFvPatchScalarField::nutUSpaldingWallFunctionFvPatchScalarField().
|
inlinevirtual |
Definition at line 241 of file nutUSpaldingWallFunctionFvPatchScalarField.H.
References nutUSpaldingWallFunctionFvPatchScalarField::nutUSpaldingWallFunctionFvPatchScalarField().
|
virtual |
Implements nutWallFunctionFvPatchScalarField.
Definition at line 311 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
References Foam::constant::atomic::group, IOobject::groupName(), Foam::mag(), turbulenceModel::nu(), Foam::foamVersion::patch, turbulenceModel::propertiesName, fvPatchField::snGrad(), U, y, and turbulenceModel::y().
|
virtual |
Reimplemented from nutWallFunctionFvPatchScalarField.
Definition at line 333 of file nutUSpaldingWallFunctionFvPatchScalarField.C.
References os(), and fvPatchField::write().
|
protected |
Definition at line 150 of file nutUSpaldingWallFunctionFvPatchScalarField.H.
|
protected |
Definition at line 153 of file nutUSpaldingWallFunctionFvPatchScalarField.H.
Referenced by nutUSpaldingWallFunctionFvPatchScalarField::calcNut().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.