This boundary condition provides a turbulence damping function, f, wall function condition for low- and high Reynolds number, turbulent flow cases. More...
Public Member Functions | |
TypeName ("fWallFunction") | |
Runtime type information. More... | |
fWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
fWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
fWallFunctionFvPatchScalarField (const fWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
Construct by mapping given fWallFunctionFvPatchScalarField. More... | |
fWallFunctionFvPatchScalarField (const fWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
virtual tmp< fvPatchScalarField > | clone () const |
Construct and return a clone. More... | |
fWallFunctionFvPatchScalarField (const fWallFunctionFvPatchScalarField &, 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 void | updateCoeffs () |
Update the coefficients associated with the patch field. More... | |
virtual void | evaluate (const Pstream::commsTypes) |
Evaluate the patchField. More... | |
virtual void | write (Ostream &) const |
Write. More... | |
![]() | |
TypeName ("fixedValue") | |
Runtime type information. More... | |
fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
Construct by mapping the given fixedValueFvPatchField<Type> More... | |
fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &) | |
Construct as copy. More... | |
fixedValueFvPatchField (const fixedValueFvPatchField< scalar > &, const DimensionedField< scalar, volMesh > &) | |
Construct as copy setting internal field reference. More... | |
virtual tmp< fvPatchField< scalar > > | clone () const |
Construct and return a clone. More... | |
virtual bool | fixesValue () const |
Return true if this patch field fixes a value. More... | |
virtual tmp< Field< scalar > > | valueInternalCoeffs (const tmp< scalarField > &) const |
Return the matrix diagonal coefficients corresponding to the. More... | |
virtual tmp< Field< scalar > > | valueBoundaryCoeffs (const tmp< scalarField > &) const |
Return the matrix source coefficients corresponding to the. More... | |
virtual tmp< Field< scalar > > | gradientInternalCoeffs () const |
Return the matrix diagonal coefficients corresponding to the. More... | |
virtual tmp< Field< scalar > > | gradientBoundaryCoeffs () const |
Return the matrix source coefficients corresponding to the. More... | |
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") | |
Runtime type information. More... | |
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 > &) | |
Construct from patch and internal field. More... | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const word &patchType) | |
Construct from patch and internal field and patch type. More... | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const Field< Type > &) | |
Construct from patch and internal field and patch field. More... | |
fvPatchField (const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &, const bool valueRequired=false) | |
Construct from patch, internal field and dictionary. More... | |
fvPatchField (const fvPatchField< Type > &, const fvPatch &, const DimensionedField< Type, volMesh > &, const fvPatchFieldMapper &) | |
Construct by mapping the given fvPatchField onto a new patch. More... | |
fvPatchField (const fvPatchField< Type > &) | |
Construct as copy. More... | |
fvPatchField (const fvPatchField< Type > &, const DimensionedField< Type, volMesh > &) | |
Construct as copy setting internal field reference. More... | |
virtual tmp< fvPatchField< Type > > | clone (const DimensionedField< Type, volMesh > &iF) const |
Construct and return a clone setting internal field reference. More... | |
Foam::tmp< Foam::fvPatchField< Type > > | NewCalculatedType (const fvPatch &p) |
Foam::tmp< Foam::fvPatchField< Type > > | NewCalculatedType (const fvPatchField< Type2 > &pf) |
virtual | ~fvPatchField () |
Destructor. More... | |
const objectRegistry & | db () const |
Return local objectRegistry. More... | |
const fvPatch & | patch () const |
Return patch. More... | |
const DimensionedField< Type, volMesh > & | dimensionedInternalField () const |
Return dimensioned internal field reference. More... | |
const Field< Type > & | internalField () const |
Return internal field reference. More... | |
const word & | patchType () const |
Optional patch type. More... | |
word & | patchType () |
Optional patch type. More... | |
virtual bool | coupled () const |
Return true if this patch field is coupled. More... | |
bool | updated () const |
Return true if the boundary condition has already been updated. More... | |
bool | manipulatedMatrix () const |
Return true if the matrix has already been manipulated. More... | |
virtual void | autoMap (const fvPatchFieldMapper &) |
Map (and resize as needed) from self given a mapping object. More... | |
virtual void | rmap (const fvPatchField< Type > &, const labelList &) |
Reverse map the given fvPatchField onto this fvPatchField. More... | |
virtual tmp< Field< Type > > | snGrad () const |
Return patch-normal gradient. More... | |
virtual tmp< Field< Type > > | snGrad (const scalarField &deltaCoeffs) const |
Return patch-normal gradient for coupled-patches. More... | |
virtual void | updateCoeffs (const scalarField &weights) |
Update the coefficients associated with the patch field. More... | |
virtual tmp< Field< Type > > | patchInternalField () const |
Return internal field next to patch as patch field. More... | |
virtual void | patchInternalField (Field< Type > &) const |
Return internal field next to patch as patch field. More... | |
virtual tmp< Field< Type > > | patchNeighbourField () const |
Return patchField on the opposite patch of a coupled patch. More... | |
virtual void | initEvaluate (const Pstream::commsTypes commsType=Pstream::blocking) |
Initialise the evaluation of the patch field. More... | |
virtual tmp< Field< Type > > | gradientInternalCoeffs (const scalarField &deltaCoeffs) const |
Return the matrix diagonal coefficients corresponding to the. More... | |
virtual tmp< Field< Type > > | gradientBoundaryCoeffs (const scalarField &deltaCoeffs) const |
Return the matrix source coefficients corresponding to the. More... | |
virtual void | manipulateMatrix (fvMatrix< Type > &matrix) |
Manipulate matrix. More... | |
virtual void | manipulateMatrix (fvMatrix< Type > &matrix, const scalarField &weights) |
Manipulate matrix with given weights. More... | |
void | writeEntryIfDifferent (Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2) const |
Helper function to write the keyword and entry only if the. More... | |
void | check (const fvPatchField< Type > &) const |
Check fvPatchField<Type> against given fvPatchField<Type> More... | |
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 | checkType () |
Check the type of the patch. More... | |
virtual void | writeLocalEntries (Ostream &) const |
Write local wall function variables. More... | |
scalar | yPlusLam (const scalar kappa, const scalar E) |
Calculate the Y+ at the edge of the laminar sublayer. More... | |
Protected Attributes | |
scalar | Cmu_ |
Cmu coefficient. More... | |
scalar | kappa_ |
Von Karman constant. More... | |
scalar | E_ |
E coefficient. More... | |
scalar | yPlusLam_ |
Y+ at the edge of the laminar sublayer. More... | |
Additional Inherited Members | |
![]() | |
typedef fvPatch | Patch |
![]() | |
static tmp< fvPatchField< Type > > | New (const word &, const fvPatch &, const DimensionedField< Type, volMesh > &) |
Return a pointer to a new patchField created on freestore given. More... | |
static tmp< fvPatchField< Type > > | New (const word &, const word &actualPatchType, const fvPatch &, const DimensionedField< Type, volMesh > &) |
Return a pointer to a new patchField created on freestore given. More... | |
static tmp< fvPatchField< Type > > | New (const fvPatchField< Type > &, const fvPatch &, const DimensionedField< Type, volMesh > &, const fvPatchFieldMapper &) |
Return a pointer to a new patchField created on freestore from. More... | |
static tmp< fvPatchField< Type > > | New (const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &) |
Return a pointer to a new patchField created on freestore. More... | |
static tmp< fvPatchField< Type > > | NewCalculatedType (const fvPatch &) |
Return a pointer to a new calculatedFvPatchField created on. More... | |
static tmp< fvPatchField< Type > > | NewCalculatedType (const fvPatchField< Type2 > &) |
Return a pointer to a new calculatedFvPatchField created on. More... | |
static const word & | calculatedType () |
Return the type of the calculated for of fvPatchField. More... | |
![]() | |
static int | disallowGenericFvPatchField |
Debug switch to disallow the use of genericFvPatchField. More... | |
This boundary condition provides a turbulence damping function, f, wall function condition for low- and high Reynolds number, turbulent flow cases.
The model operates in two modes, based on the computed laminar-to-turbulent switch-over y+ value derived from kappa and E.
Patch usage
Property | Description | Required | Default value |
---|---|---|---|
Cmu | model coefficient | no | 0.09 |
kappa | Von Karman constant | no | 0.41 |
E | model coefficient | no | 9.8 |
Example of the boundary condition specification:
myPatch { type fWallFunction; }
Definition at line 97 of file fWallFunctionFvPatchScalarField.H.
fWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 84 of file fWallFunctionFvPatchScalarField.C.
Referenced by fWallFunctionFvPatchScalarField::clone().
fWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Construct from patch, internal field and dictionary.
Definition at line 118 of file fWallFunctionFvPatchScalarField.C.
fWallFunctionFvPatchScalarField | ( | const fWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Construct by mapping given fWallFunctionFvPatchScalarField.
onto a new patch
Definition at line 100 of file fWallFunctionFvPatchScalarField.C.
fWallFunctionFvPatchScalarField | ( | const fWallFunctionFvPatchScalarField & | v2wfpsf | ) |
Construct as copy.
Definition at line 135 of file fWallFunctionFvPatchScalarField.C.
fWallFunctionFvPatchScalarField | ( | const fWallFunctionFvPatchScalarField & | v2wfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct as copy setting internal field reference.
Definition at line 150 of file fWallFunctionFvPatchScalarField.C.
|
protectedvirtual |
Check the type of the patch.
Definition at line 42 of file fWallFunctionFvPatchScalarField.C.
References Foam::abort(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, fvPatch::name(), Foam::nl, and fvPatchField< scalar >::patch().
|
protectedvirtual |
Write local wall function variables.
Definition at line 56 of file fWallFunctionFvPatchScalarField.C.
References fWallFunctionFvPatchScalarField::Cmu_, fWallFunctionFvPatchScalarField::E_, token::END_STATEMENT, fWallFunctionFvPatchScalarField::kappa_, Foam::nl, and Ostream::writeKeyword().
Referenced by fWallFunctionFvPatchScalarField::write().
|
protected |
Calculate the Y+ at the edge of the laminar sublayer.
Definition at line 65 of file fWallFunctionFvPatchScalarField.C.
References Foam::constant::electromagnetic::kappa, Foam::log(), and Foam::max().
TypeName | ( | "fWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Construct and return a clone.
Definition at line 170 of file fWallFunctionFvPatchScalarField.H.
References fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField().
|
inlinevirtual |
Construct and return a clone setting internal field reference.
Reimplemented from fixedValueFvPatchField< scalar >.
Definition at line 187 of file fWallFunctionFvPatchScalarField.H.
References fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField().
|
virtual |
Update the coefficients associated with the patch field.
Reimplemented from fvPatchField< scalar >.
Definition at line 167 of file fWallFunctionFvPatchScalarField.C.
References fWallFunctionFvPatchScalarField::Cmu_, fvPatchField< scalar >::db(), fvPatchField< scalar >::dimensionedInternalField(), epsilon, turbulenceModel::epsilon(), f(), fvPatch::faceCells(), forAll, Foam::constant::atomic::group, IOobject::groupName(), fvPatch::index(), k, turbulenceModel::k(), objectRegistry::lookupObject(), N, turbulenceModel::nu(), fvPatchField< scalar >::patch(), patchi, Foam::pow025(), turbulenceModel::propertiesName, Foam::sqr(), Foam::sqrt(), fvPatchField< Type >::updateCoeffs(), fvPatchField< scalar >::updated(), uTau, y, turbulenceModel::y(), and fWallFunctionFvPatchScalarField::yPlusLam_.
|
virtual |
Evaluate the patchField.
Reimplemented from fvPatchField< scalar >.
Definition at line 236 of file fWallFunctionFvPatchScalarField.C.
References fvPatchField< Type >::evaluate().
|
virtual |
Write.
Reimplemented from fixedValueFvPatchField< scalar >.
Definition at line 244 of file fWallFunctionFvPatchScalarField.C.
References fixedValueFvPatchField< Type >::write(), and fWallFunctionFvPatchScalarField::writeLocalEntries().
|
protected |
Cmu coefficient.
Definition at line 106 of file fWallFunctionFvPatchScalarField.H.
Referenced by fWallFunctionFvPatchScalarField::updateCoeffs(), and fWallFunctionFvPatchScalarField::writeLocalEntries().
|
protected |
Von Karman constant.
Definition at line 109 of file fWallFunctionFvPatchScalarField.H.
Referenced by fWallFunctionFvPatchScalarField::writeLocalEntries().
|
protected |
E coefficient.
Definition at line 112 of file fWallFunctionFvPatchScalarField.H.
Referenced by fWallFunctionFvPatchScalarField::writeLocalEntries().
|
protected |
Y+ at the edge of the laminar sublayer.
Definition at line 115 of file fWallFunctionFvPatchScalarField.H.
Referenced by fWallFunctionFvPatchScalarField::updateCoeffs().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.