This boundary condition provides a turbulence dissipation wall function condition for high Reynolds number, turbulent flow cases. More...
Public Member Functions | |
TypeName ("epsilonWallFunction") | |
Runtime type information. More... | |
epsilonWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
epsilonWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
epsilonWallFunctionFvPatchScalarField (const epsilonWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
Construct by mapping given. More... | |
epsilonWallFunctionFvPatchScalarField (const epsilonWallFunctionFvPatchScalarField &) | |
Construct as copy. More... | |
virtual tmp< fvPatchScalarField > | clone () const |
Construct and return a clone. More... | |
epsilonWallFunctionFvPatchScalarField (const epsilonWallFunctionFvPatchScalarField &, 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 | ~epsilonWallFunctionFvPatchScalarField () |
Destructor. More... | |
scalarField & | G (bool init=false) |
Return non-const access to the master's G field. More... | |
scalarField & | epsilon (bool init=false) |
Return non-const access to the master's epsilon field. More... | |
virtual void | updateCoeffs () |
Update the coefficients associated with the patch field. More... | |
virtual void | updateCoeffs (const scalarField &weights) |
Update the coefficients associated with the patch field. More... | |
virtual void | manipulateMatrix (fvMatrix< scalar > &matrix) |
Manipulate matrix. More... | |
virtual void | manipulateMatrix (fvMatrix< scalar > &matrix, const scalarField &weights) |
Manipulate matrix with given weights. 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 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 void | evaluate (const Pstream::commsTypes commsType=Pstream::blocking) |
Evaluate the patch field, sets Updated to false. 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... | |
virtual void | setMaster () |
Set the master patch - master is responsible for updating all. More... | |
virtual void | createAveragingWeights () |
Create the averaging weights for cells which are bounded by. More... | |
virtual epsilonWallFunctionFvPatchScalarField & | epsilonPatch (const label patchi) |
Helper function to return non-const access to an epsilon patch. More... | |
virtual void | calculateTurbulenceFields (const turbulenceModel &turbulence, scalarField &G0, scalarField &epsilon0) |
Main driver to calculate the turbulence fields. More... | |
virtual void | calculate (const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon) |
Calculate the epsilon and G. More... | |
virtual label & | master () |
Return non-const access to the master patch ID. More... | |
Protected Attributes | |
scalar | Cmu_ |
Cmu coefficient. More... | |
scalar | kappa_ |
Von Karman constant. More... | |
scalar | E_ |
E coefficient. More... | |
scalarField | G_ |
Local copy of turbulence G field. More... | |
scalarField | epsilon_ |
Local copy of turbulence epsilon field. More... | |
bool | initialised_ |
Initialised flag. More... | |
label | master_ |
Master patch ID. More... | |
List< List< scalar > > | cornerWeights_ |
List of averaging corner weights. More... | |
Static Protected Attributes | |
static scalar | tolerance_ = 1e-5 |
Tolerance used in weighted calculations. 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 dissipation wall function condition for high Reynolds number, turbulent flow cases.
The condition can be applied to wall boundaries, whereby it
epsilon
and G
where
![]() | = | turblence dissipation field |
![]() | = | turblence generation field |
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 epsilonWallFunction; }
Definition at line 114 of file epsilonWallFunctionFvPatchScalarField.H.
epsilonWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 254 of file epsilonWallFunctionFvPatchScalarField.C.
Referenced by epsilonWallFunctionFvPatchScalarField::clone().
epsilonWallFunctionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Construct from patch, internal field and dictionary.
Definition at line 298 of file epsilonWallFunctionFvPatchScalarField.C.
References Foam::operator==().
epsilonWallFunctionFvPatchScalarField | ( | const epsilonWallFunctionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Construct by mapping given.
epsilonWallFunctionFvPatchScalarField onto a new patch
Definition at line 275 of file epsilonWallFunctionFvPatchScalarField.C.
epsilonWallFunctionFvPatchScalarField | ( | const epsilonWallFunctionFvPatchScalarField & | ewfpsf | ) |
Construct as copy.
Definition at line 323 of file epsilonWallFunctionFvPatchScalarField.C.
epsilonWallFunctionFvPatchScalarField | ( | const epsilonWallFunctionFvPatchScalarField & | ewfpsf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct as copy setting internal field reference.
Definition at line 343 of file epsilonWallFunctionFvPatchScalarField.C.
|
inlinevirtual |
Destructor.
Definition at line 266 of file epsilonWallFunctionFvPatchScalarField.H.
|
protectedvirtual |
Check the type of the patch.
Definition at line 40 of file epsilonWallFunctionFvPatchScalarField.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 55 of file epsilonWallFunctionFvPatchScalarField.C.
References token::END_STATEMENT, Foam::nl, and Ostream::writeKeyword().
|
protectedvirtual |
Set the master patch - master is responsible for updating all.
wall function patches
Definition at line 65 of file epsilonWallFunctionFvPatchScalarField.C.
References dimensionedInternalField(), epsilon, forAll, epsilonWallFunctionFvPatchScalarField::master(), and patchi.
|
protectedvirtual |
Create the averaging weights for cells which are bounded by.
multiple wall function faces
Definition at line 95 of file epsilonWallFunctionFvPatchScalarField.C.
References DynamicList::append(), GeometricField::boundaryField(), polyMesh::changing(), dimensionedInternalField(), Foam::dimless, epsilon, forAll, mesh, IOobject::NO_READ, IOobject::NO_WRITE, patchi, fvPatchField::patchInternalField(), fvMesh::time(), and Time::timeName().
|
protectedvirtual |
Helper function to return non-const access to an epsilon patch.
Definition at line 155 of file epsilonWallFunctionFvPatchScalarField.C.
References dimensionedInternalField(), epsilon, and patchi.
|
protectedvirtual |
Main driver to calculate the turbulence fields.
Definition at line 170 of file epsilonWallFunctionFvPatchScalarField.C.
References epsilonWallFunctionFvPatchScalarField::calculate(), Foam::constant::electromagnetic::epsilon0, fvPatch::faceCells(), forAll, Foam::constant::electromagnetic::G0, fvPatchField::patch(), patchi, turbulence, and w().
|
protectedvirtual |
Calculate the epsilon and G.
Reimplemented in epsilonLowReWallFunctionFvPatchScalarField.
Definition at line 203 of file epsilonWallFunctionFvPatchScalarField.C.
References epsilon, fvPatch::faceCells(), forAll, Foam::constant::universal::G, fvPatch::index(), k, Foam::mag(), patchi, Foam::pow(), Foam::pow025(), fvPatchField::snGrad(), Foam::sqrt(), turbulence, w(), and y.
Referenced by epsilonWallFunctionFvPatchScalarField::calculateTurbulenceFields().
|
inlineprotectedvirtual |
Return non-const access to the master patch ID.
Definition at line 191 of file epsilonWallFunctionFvPatchScalarField.H.
References epsilonWallFunctionFvPatchScalarField::master_.
Referenced by epsilonWallFunctionFvPatchScalarField::setMaster().
TypeName | ( | "epsilonWallFunction" | ) |
Runtime type information.
|
inlinevirtual |
Construct and return a clone.
Reimplemented in epsilonLowReWallFunctionFvPatchScalarField.
Definition at line 238 of file epsilonWallFunctionFvPatchScalarField.H.
References epsilonWallFunctionFvPatchScalarField::epsilonWallFunctionFvPatchScalarField().
|
inlinevirtual |
Construct and return a clone setting internal field reference.
Reimplemented from fixedValueFvPatchField< scalar >.
Reimplemented in epsilonLowReWallFunctionFvPatchScalarField.
Definition at line 255 of file epsilonWallFunctionFvPatchScalarField.H.
References epsilonWallFunctionFvPatchScalarField::epsilonWallFunctionFvPatchScalarField().
Foam::scalarField & G | ( | bool | init = false | ) |
Return non-const access to the master's G field.
Definition at line 364 of file epsilonWallFunctionFvPatchScalarField.C.
Foam::scalarField & epsilon | ( | bool | init = false | ) |
Return non-const access to the master's epsilon field.
Definition at line 381 of file epsilonWallFunctionFvPatchScalarField.C.
|
virtual |
Update the coefficients associated with the patch field.
Reimplemented from fvPatchField< scalar >.
Definition at line 399 of file epsilonWallFunctionFvPatchScalarField.C.
References dimensionedInternalField(), epsilon, Foam::constant::electromagnetic::epsilon0, fvPatch::faceCells(), forAll, Foam::constant::universal::G, Foam::constant::electromagnetic::G0, turbulenceModel::GName(), Foam::constant::atomic::group, IOobject::groupName(), turbulenceModel::propertiesName, and fvPatchField::updateCoeffs().
|
virtual |
Update the coefficients associated with the patch field.
Reimplemented from fvPatchField< scalar >.
Definition at line 449 of file epsilonWallFunctionFvPatchScalarField.C.
References dimensionedInternalField(), epsilon, Foam::constant::electromagnetic::epsilon0, fvPatch::faceCells(), forAll, Foam::constant::universal::G, Foam::constant::electromagnetic::G0, turbulenceModel::GName(), Foam::constant::atomic::group, IOobject::groupName(), turbulenceModel::propertiesName, fvPatchField::updateCoeffs(), and w().
|
virtual |
Manipulate matrix.
Definition at line 510 of file epsilonWallFunctionFvPatchScalarField.C.
References fvPatchField::manipulateMatrix(), and fvMatrix::setValues().
|
virtual |
Manipulate matrix with given weights.
Definition at line 526 of file epsilonWallFunctionFvPatchScalarField.C.
References DynamicList::append(), dimensionedInternalField(), Foam::endl(), epsilon, fvPatch::faceCells(), forAll, fvPatchField::manipulateMatrix(), fvPatch::name(), Foam::Pout, fvMatrix::setValues(), fvPatch::size(), and DynamicList::xfer().
|
virtual |
Write.
Reimplemented from fixedValueFvPatchField< scalar >.
Definition at line 578 of file epsilonWallFunctionFvPatchScalarField.C.
References fixedValueFvPatchField< Type >::write().
|
staticprotected |
Tolerance used in weighted calculations.
Definition at line 123 of file epsilonWallFunctionFvPatchScalarField.H.
|
protected |
Cmu coefficient.
Definition at line 126 of file epsilonWallFunctionFvPatchScalarField.H.
|
protected |
Von Karman constant.
Definition at line 129 of file epsilonWallFunctionFvPatchScalarField.H.
|
protected |
E coefficient.
Definition at line 132 of file epsilonWallFunctionFvPatchScalarField.H.
|
protected |
Local copy of turbulence G field.
Definition at line 135 of file epsilonWallFunctionFvPatchScalarField.H.
|
protected |
Local copy of turbulence epsilon field.
Definition at line 138 of file epsilonWallFunctionFvPatchScalarField.H.
|
protected |
Initialised flag.
Definition at line 141 of file epsilonWallFunctionFvPatchScalarField.H.
|
protected |
Master patch ID.
Definition at line 144 of file epsilonWallFunctionFvPatchScalarField.H.
Referenced by epsilonWallFunctionFvPatchScalarField::master().
List of averaging corner weights.
Definition at line 147 of file epsilonWallFunctionFvPatchScalarField.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.