This boundary condition provides a wall constraint on the specific dissipation rate (i.e. omega
) and the turbulent kinetic energy production contribution (i.e. G
) for atmospheric boundary layer modelling.
More...
Public Member Functions | |
TypeName ("atmOmegaWallFunction") | |
atmOmegaWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
atmOmegaWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
atmOmegaWallFunctionFvPatchScalarField (const atmOmegaWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
atmOmegaWallFunctionFvPatchScalarField (const atmOmegaWallFunctionFvPatchScalarField &) | |
virtual tmp< fvPatchScalarField > | clone () const |
atmOmegaWallFunctionFvPatchScalarField (const atmOmegaWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
virtual tmp< fvPatchScalarField > | clone (const DimensionedField< scalar, volMesh > &iF) const |
virtual | ~atmOmegaWallFunctionFvPatchScalarField ()=default |
virtual void | autoMap (const fvPatchFieldMapper &) |
virtual void | rmap (const fvPatchScalarField &, const labelList &) |
virtual void | write (Ostream &) const |
![]() | |
TypeName ("omegaWallFunction") | |
omegaWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
omegaWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
omegaWallFunctionFvPatchScalarField (const omegaWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
omegaWallFunctionFvPatchScalarField (const omegaWallFunctionFvPatchScalarField &) | |
omegaWallFunctionFvPatchScalarField (const omegaWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
virtual | ~omegaWallFunctionFvPatchScalarField ()=default |
scalarField & | G (bool init=false) |
scalarField & | omega (bool init=false) |
virtual void | updateCoeffs () |
virtual void | updateWeightedCoeffs (const scalarField &weights) |
virtual void | manipulateMatrix (fvMatrix< scalar > &matrix) |
virtual void | manipulateMatrix (fvMatrix< scalar > &matrix, const scalarField &weights) |
![]() | |
TypeName ("fixedValue") | |
fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const scalar &value) | |
fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &, const bool valueRequired=true) | |
fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &) | |
fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &, const DimensionedField< scalar, volMesh > &) | |
virtual tmp< fvPatchField< scalar > > | clone () const |
virtual bool | fixesValue () const |
virtual bool | assignable () const |
virtual tmp< Field< scalar > > | valueInternalCoeffs (const tmp< scalarField > &) const |
virtual tmp< Field< scalar > > | valueBoundaryCoeffs (const tmp< scalarField > &) const |
virtual tmp< Field< scalar > > | gradientInternalCoeffs () const |
virtual tmp< Field< scalar > > | gradientBoundaryCoeffs () const |
virtual void | operator= (const UList< scalar > &) |
virtual void | operator= (const fvPatchField< scalar > &) |
virtual void | operator= (const scalar &) |
virtual void | operator+= (const fvPatchField< scalar > &) |
virtual void | operator+= (const Field< scalar > &) |
virtual void | operator+= (const scalar &) |
virtual void | operator-= (const fvPatchField< scalar > &) |
virtual void | operator-= (const Field< scalar > &) |
virtual void | operator-= (const scalar &) |
virtual void | operator*= (const fvPatchField< scalar > &) |
virtual void | operator*= (const Field< scalar > &) |
virtual void | operator*= (const scalar) |
virtual void | operator/= (const fvPatchField< scalar > &) |
virtual void | operator/= (const Field< scalar > &) |
virtual void | operator/= (const scalar) |
![]() | |
TypeName ("fvPatchField") | |
declareRunTimeSelectionTable (tmp, fvPatchField, patch,(const fvPatch &p, const DimensionedField< Type, volMesh > &iF),(p, iF)) | |
declareRunTimeSelectionTable (tmp, fvPatchField, patchMapper,(const fvPatchField< Type > &ptf, const fvPatch &p, const DimensionedField< Type, volMesh > &iF, const fvPatchFieldMapper &m),(dynamic_cast< const fvPatchFieldType & >(ptf), p, iF, m)) | |
declareRunTimeSelectionTable (tmp, fvPatchField, dictionary,(const fvPatch &p, const DimensionedField< Type, volMesh > &iF, const dictionary &dict),(p, iF, dict)) | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &) | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const Type &value) | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const word &patchType) | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const Field< Type > &) | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &, const bool valueRequired=true) | |
fvPatchField (const fvPatchField< Type > &, const fvPatch &, const DimensionedField< Type, volMesh > &, const fvPatchFieldMapper &) | |
fvPatchField (const fvPatchField< Type > &) | |
fvPatchField (const fvPatchField< Type > &, const DimensionedField< Type, volMesh > &) | |
virtual tmp< fvPatchField< Type > > | clone (const DimensionedField< Type, volMesh > &iF) const |
Foam::tmp< Foam::fvPatchField< Type > > | NewCalculatedType (const fvPatch &p) |
Foam::tmp< Foam::fvPatchField< Type > > | NewCalculatedType (const fvPatchField< Type2 > &pf) |
virtual | ~fvPatchField ()=default |
bool | useImplicit () const noexcept |
bool | useImplicit (bool on) noexcept |
virtual bool | coupled () const |
const objectRegistry & | db () const |
const fvPatch & | patch () const |
const DimensionedField< Type, volMesh > & | internalField () const |
const Field< Type > & | primitiveField () const |
const word & | patchType () const |
word & | patchType () |
bool | updated () const |
bool | manipulatedMatrix () const |
virtual void | rmap (const fvPatchField< Type > &, const labelList &) |
virtual tmp< Field< Type > > | snGrad () const |
virtual tmp< Field< Type > > | snGrad (const scalarField &deltaCoeffs) const |
virtual tmp< Field< Type > > | patchInternalField () const |
virtual void | patchInternalField (Field< Type > &) const |
virtual tmp< Field< Type > > | patchNeighbourField () const |
virtual void | initEvaluate (const Pstream::commsTypes commsType=Pstream::commsTypes::blocking) |
virtual void | evaluate (const Pstream::commsTypes commsType=Pstream::commsTypes::blocking) |
virtual tmp< Field< Type > > | valueInternalCoeffs (const tmp< Field< scalar >> &) const |
virtual tmp< Field< Type > > | valueBoundaryCoeffs (const tmp< Field< scalar >> &) const |
virtual tmp< Field< Type > > | gradientInternalCoeffs (const scalarField &deltaCoeffs) const |
virtual tmp< Field< Type > > | gradientBoundaryCoeffs (const scalarField &deltaCoeffs) const |
virtual void | manipulateMatrix (fvMatrix< Type > &matrix) |
virtual void | manipulateMatrix (fvMatrix< Type > &matrix, const scalarField &weights) |
virtual void | manipulateMatrix (fvMatrix< Type > &matrix, const label iMatrix, const direction cmp) |
void | check (const fvPatchField< Type > &) const |
virtual void | operator= (const UList< Type > &) |
virtual void | operator= (const fvPatchField< Type > &) |
virtual void | operator= (const Type &) |
virtual void | operator+= (const fvPatchField< Type > &) |
virtual void | operator+= (const Field< Type > &) |
virtual void | operator+= (const Type &) |
virtual void | operator-= (const fvPatchField< Type > &) |
virtual void | operator-= (const Field< Type > &) |
virtual void | operator-= (const Type &) |
virtual void | operator== (const fvPatchField< Type > &) |
virtual void | operator== (const Field< Type > &) |
virtual void | operator== (const Type &) |
Protected Member Functions | |
virtual void | calculate (const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega) |
![]() | |
virtual void | setMaster () |
virtual void | createAveragingWeights () |
virtual omegaWallFunctionFvPatchScalarField & | omegaPatch (const label patchi) |
virtual void | calculateTurbulenceFields (const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0) |
virtual label & | master () |
Protected Attributes | |
autoPtr< PatchFunction1< scalar > > | z0_ |
![]() | |
bool | blended_ |
bool | initialised_ |
label | master_ |
scalar | beta1_ |
scalarField | G_ |
scalarField | omega_ |
List< List< scalar > > | cornerWeights_ |
Additional Inherited Members | |
![]() | |
typedef fvPatch | Patch |
typedef calculatedFvPatchField< Type > | Calculated |
![]() | |
static tmp< fvPatchField< Type > > | New (const word &, const fvPatch &, const DimensionedField< Type, volMesh > &) |
static tmp< fvPatchField< Type > > | New (const word &, const word &actualPatchType, const fvPatch &, const DimensionedField< Type, volMesh > &) |
static tmp< fvPatchField< Type > > | New (const fvPatchField< Type > &, const fvPatch &, const DimensionedField< Type, volMesh > &, const fvPatchFieldMapper &) |
static tmp< fvPatchField< Type > > | New (const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &) |
static tmp< fvPatchField< Type > > | NewCalculatedType (const fvPatch &) |
static tmp< fvPatchField< Type > > | NewCalculatedType (const fvPatchField< Type2 > &) |
static const word & | calculatedType () |
![]() | |
static int | disallowGenericFvPatchField |
![]() | |
static scalar | tolerance_ = 1e-5 |
This boundary condition provides a wall constraint on the specific dissipation rate (i.e. omega
) and the turbulent kinetic energy production contribution (i.e. G
) for atmospheric boundary layer modelling.
Theoretical expressions (tags:PGVB, B): Parente, A., Gorlé, C., Van Beeck, J., & Benocci, C. (2011). Improved k–ε model and wall function formulation for the RANS simulation of ABL flows. J. of wind engineering and industrial aerodynamics, 99(4), 267-278. DOI:10.1016/j.jweia.2010.12.017 Bredberg, J. (2000). On the wall boundary condition for turbulence models. Chalmers University of Technology, Depart. of Thermo and Fluid Dyn. Internal Report 00/4. Sweden: Göteborg.
Required fields:
omega | Specific dissipation rate [1/s]
<patchName> { // Mandatory entries (unmodifiable) type atmOmegaWallFunction; // Mandatory entries (runtime modifiable) z0 uniform 0.001; // Optional (inherited) entries ... }
where the entries mean:
Property | Description | Type | Reqd | Dflt |
---|---|---|---|---|
type | Type name: atmOmegaWallFunction | word | yes | - |
z0 | Surface roughness length [m] | PatchFunction1<scalar> | yes | - |
The inherited entries are elaborated in:
Definition at line 122 of file atmOmegaWallFunctionFvPatchScalarField.H.
atmOmegaWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Definition at line 100 of file atmOmegaWallFunctionFvPatchScalarField.C.
Referenced by atmOmegaWallFunctionFvPatchScalarField::clone().
atmOmegaWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Definition at line 126 of file atmOmegaWallFunctionFvPatchScalarField.C.
atmOmegaWallFunctionFvPatchScalarField | ( | const atmOmegaWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Definition at line 112 of file atmOmegaWallFunctionFvPatchScalarField.C.
atmOmegaWallFunctionFvPatchScalarField | ( | const atmOmegaWallFunctionFvPatchScalarField & | owfpsf | ) |
Definition at line 139 of file atmOmegaWallFunctionFvPatchScalarField.C.
atmOmegaWallFunctionFvPatchScalarField | ( | const atmOmegaWallFunctionFvPatchScalarField & | owfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Definition at line 150 of file atmOmegaWallFunctionFvPatchScalarField.C.
|
virtualdefault |
|
protectedvirtual |
Reimplemented from omegaWallFunctionFvPatchScalarField.
Definition at line 31 of file atmOmegaWallFunctionFvPatchScalarField.C.
References GeometricField::boundaryField(), nutWallFunctionFvPatchScalarField::Cmu(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::constant::electromagnetic::G0, k, turbulenceModel::k(), nutWallFunctionFvPatchScalarField::kappa(), Foam::mag(), turbulenceModel::nu(), Foam::foamVersion::patch, Foam::pow025(), fvPatchField::snGrad(), Foam::sqrt(), turbulenceModel::U(), y, and turbulenceModel::y().
TypeName | ( | "atmOmegaWallFunction" | ) |
|
inlinevirtual |
Reimplemented from omegaWallFunctionFvPatchScalarField.
Definition at line 188 of file atmOmegaWallFunctionFvPatchScalarField.H.
References atmOmegaWallFunctionFvPatchScalarField::atmOmegaWallFunctionFvPatchScalarField().
|
inlinevirtual |
Reimplemented from omegaWallFunctionFvPatchScalarField.
Definition at line 205 of file atmOmegaWallFunctionFvPatchScalarField.H.
References atmOmegaWallFunctionFvPatchScalarField::atmOmegaWallFunctionFvPatchScalarField().
|
virtual |
Reimplemented from fvPatchField< scalar >.
Definition at line 163 of file atmOmegaWallFunctionFvPatchScalarField.C.
|
virtual |
Definition at line 173 of file atmOmegaWallFunctionFvPatchScalarField.C.
References atmOmegaWallFunctionFvPatchScalarField::z0_.
|
virtual |
Reimplemented from omegaWallFunctionFvPatchScalarField.
Definition at line 188 of file atmOmegaWallFunctionFvPatchScalarField.C.
References os(), and Foam::vtk::write().
|
protected |
Definition at line 131 of file atmOmegaWallFunctionFvPatchScalarField.H.
Referenced by atmOmegaWallFunctionFvPatchScalarField::rmap().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.