omegaWallFunctionFvPatchScalarField Class Reference

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 low- and high-Reynolds number turbulence models. More...

Inheritance diagram for omegaWallFunctionFvPatchScalarField:
[legend]
Collaboration diagram for omegaWallFunctionFvPatchScalarField:
[legend]

Public Member Functions

 TypeName ("omegaWallFunction")
 Runtime type information. More...
 
 omegaWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Construct from patch and internal field. More...
 
 omegaWallFunctionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 Construct from patch, internal field and dictionary. More...
 
 omegaWallFunctionFvPatchScalarField (const omegaWallFunctionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 
 omegaWallFunctionFvPatchScalarField (const omegaWallFunctionFvPatchScalarField &)
 Construct as copy. More...
 
virtual tmp< fvPatchScalarFieldclone () const
 Construct and return a clone. More...
 
 omegaWallFunctionFvPatchScalarField (const omegaWallFunctionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &)
 Construct as copy setting internal field reference. More...
 
virtual tmp< fvPatchScalarFieldclone (const DimensionedField< scalar, volMesh > &iF) const
 Construct and return a clone setting internal field reference. More...
 
virtual ~omegaWallFunctionFvPatchScalarField ()=default
 Destructor. More...
 
scalarFieldG (bool init=false)
 Return non-const access to the master's G field. More...
 
scalarFieldomega (bool init=false)
 Return non-const access to the master's omega field. More...
 
virtual void updateCoeffs ()
 Update the coefficients associated with the patch field. More...
 
virtual void updateWeightedCoeffs (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...
 
- Public Member Functions inherited from fixedValueFvPatchField< scalar >
 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 scalar &value)
 Construct from patch, internal field and value. More...
 
 fixedValueFvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &, const bool valueRequired=true)
 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 bool fixesValue () const
 Return true if this patch field fixes a value. More...
 
virtual bool assignable () const
 Return false: this patch field is not altered by assignment. 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)
 
- Public Member Functions inherited from fvPatchField< scalar >
 TypeName ("fvPatchField")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (tmp, fvPatchField, patch,(const fvPatch &p, const DimensionedField< scalar, volMesh > &iF),(p, iF))
 
 declareRunTimeSelectionTable (tmp, fvPatchField, patchMapper,(const fvPatchField< scalar > &ptf, const fvPatch &p, const DimensionedField< scalar, volMesh > &iF, const fvPatchFieldMapper &m),(dynamic_cast< const fvPatchFieldType & >(ptf), p, iF, m))
 
 declareRunTimeSelectionTable (tmp, fvPatchField, dictionary,(const fvPatch &p, const DimensionedField< scalar, volMesh > &iF, const dictionary &dict),(p, iF, dict))
 
 fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Construct from patch and internal field. More...
 
 fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const scalar &value)
 Construct from patch, internal field and value. More...
 
 fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const word &patchType)
 Construct from patch and internal field and patch type. More...
 
 fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const Field< scalar > &)
 Construct from patch and internal field and patch field. More...
 
 fvPatchField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &, const bool valueRequired=true)
 Construct from patch, internal field and dictionary. More...
 
 fvPatchField (const fvPatchField< scalar > &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 Construct by mapping the given fvPatchField onto a new patch. More...
 
 fvPatchField (const fvPatchField< scalar > &)
 Construct as copy. More...
 
 fvPatchField (const fvPatchField< scalar > &, const DimensionedField< scalar, volMesh > &)
 Construct as copy setting internal field reference. More...
 
Foam::tmp< Foam::fvPatchField< scalar > > NewCalculatedType (const fvPatchField< Type2 > &pf)
 
virtual ~fvPatchField ()=default
 Destructor. More...
 
bool useImplicit () const noexcept
 Use implicit formulation for coupled patches only. More...
 
bool useImplicit (bool on) noexcept
 Set useImplicit on/off. More...
 
virtual bool coupled () const
 Return true if this patch field is coupled. More...
 
const objectRegistrydb () const
 Return local objectRegistry. More...
 
const fvPatchpatch () const
 Return patch. More...
 
const DimensionedField< scalar, volMesh > & internalField () const
 Return dimensioned internal field reference. More...
 
const Field< scalar > & primitiveField () const
 Return internal field reference. More...
 
const wordpatchType () const
 Optional patch type. More...
 
wordpatchType ()
 Optional patch type. 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< scalar > &, const labelList &)
 Reverse map the given fvPatchField onto this fvPatchField. More...
 
virtual tmp< Field< scalar > > snGrad () const
 Return patch-normal gradient. More...
 
virtual tmp< Field< scalar > > snGrad (const scalarField &deltaCoeffs) const
 Return patch-normal gradient for coupled-patches. More...
 
virtual tmp< Field< scalar > > patchInternalField () const
 Return internal field next to patch as patch field. More...
 
virtual void patchInternalField (Field< scalar > &) const
 Return internal field next to patch as patch field. More...
 
virtual tmp< Field< scalar > > patchNeighbourField () const
 Return patchField on the opposite patch of a coupled patch. More...
 
virtual void initEvaluate (const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
 Initialise the evaluation of the patch field. More...
 
virtual void evaluate (const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
 Evaluate the patch field, sets Updated to false. More...
 
virtual tmp< Field< scalar > > valueInternalCoeffs (const tmp< Field< scalar >> &) const
 Return the matrix diagonal coefficients corresponding to the. More...
 
virtual tmp< Field< scalar > > valueBoundaryCoeffs (const tmp< Field< scalar >> &) const
 Return the matrix source coefficients corresponding to the. More...
 
virtual tmp< Field< scalar > > gradientInternalCoeffs (const scalarField &deltaCoeffs) const
 Return the matrix diagonal coefficients corresponding to the. More...
 
virtual tmp< Field< scalar > > gradientBoundaryCoeffs (const scalarField &deltaCoeffs) const
 Return the matrix source coefficients corresponding to the. More...
 
virtual void manipulateMatrix (fvMatrix< scalar > &matrix, const label iMatrix, const direction cmp)
 Manipulate fvMatrix. More...
 
void check (const fvPatchField< scalar > &) const
 Check fvPatchField<Type> against given fvPatchField<Type> More...
 
virtual void operator== (const fvPatchField< scalar > &)
 
virtual void operator== (const Field< scalar > &)
 
virtual void operator== (const scalar &)
 

Protected Member Functions

virtual void setMaster ()
 
virtual void createAveragingWeights ()
 
virtual omegaWallFunctionFvPatchScalarFieldomegaPatch (const label patchi)
 Helper function to return non-const access to an omega patch. More...
 
virtual void calculateTurbulenceFields (const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
 Main driver to calculate the turbulence fields. More...
 
virtual void calculate (const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
 Calculate the omega and G. More...
 
virtual label & master ()
 Return non-const access to the master patch ID. More...
 

Protected Attributes

bool blended_
 Deprecated(2019-11) Blending switch. More...
 
bool initialised_
 Initialised flag. More...
 
label master_
 Master patch ID. More...
 
scalar beta1_
 beta1 coefficient More...
 
scalarField G_
 Local copy of turbulence G field. More...
 
scalarField omega_
 Local copy of turbulence omega field. 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

- Public Types inherited from fvPatchField< scalar >
typedef fvPatch Patch
 
typedef calculatedFvPatchField< scalar > Calculated
 
- Static Public Member Functions inherited from fvPatchField< scalar >
static tmp< fvPatchField< scalar > > New (const word &, const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Return a pointer to a new patchField created on freestore given. More...
 
static tmp< fvPatchField< scalar > > New (const word &, const word &actualPatchType, const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Return a pointer to a new patchField created on freestore given. More...
 
static tmp< fvPatchField< scalar > > New (const fvPatchField< scalar > &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 Return a pointer to a new patchField created on freestore from. More...
 
static tmp< fvPatchField< scalar > > New (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 Return a pointer to a new patchField created on freestore. More...
 
static tmp< fvPatchField< scalar > > NewCalculatedType (const fvPatch &)
 Return a pointer to a new calculatedFvPatchField created on. More...
 
static tmp< fvPatchField< scalar > > NewCalculatedType (const fvPatchField< Type2 > &)
 Return a pointer to a new calculatedFvPatchField created on. More...
 
static const wordcalculatedType ()
 Return the type of the calculated for of fvPatchField. More...
 
- Static Public Attributes inherited from fvPatchField< scalar >
static int disallowGenericFvPatchField
 Debug switch to disallow the use of genericFvPatchField. More...
 

Detailed Description

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 low- and high-Reynolds number turbulence models.

Reference:

    Binomial blending of the viscous and inertial sublayers (tag:ME):
        Menter, F., & Esch, T. (2001).
        Elements of industrial heat transfer prediction.
        In Proceedings of the 16th Brazilian Congress of Mechanical
        Engineering (COBEM), November 2001. vol. 20, p. 117-127.

    Exponential/Max blending of the viscous and inertial sublayers (tag:PH):
        Popovac, M., & Hanjalić, K. (2007).
        Compound wall treatment for RANS computation of complex
        turbulent flows and heat transfer.
        Flow, turbulence and combustion, 78(2), 177-202.
        DOI:10.1007/s10494-006-9067-x

    Tanh blending of the viscous and inertial sublayers (tag:KAS):
        Knopp, T., Alrutz, T., & Schwamborn, D. (2006).
        A grid and flow adaptive wall-function method for RANS
        turbulence modelling.
        Journal of Computational Physics, 220(1), 19-40.
        DOI:10.1016/j.jcp.2006.05.003
Usage
Example of the boundary condition specification:
<patchName>
{
    // Mandatory entries (unmodifiable)
    type            omegaWallFunction;

    // Optional entries (unmodifiable)
    beta1           0.075;
    blending        binomial2;
    n               2.0;

    // Optional (inherited) entries
    ...
}
Property Description Type Req'd Dflt
type Type name: omegaWallFunction word yes -
beta1 Model coefficient scalar no 0.075
blending Viscous/inertial sublayer blending method word no binomial2
n Binomial blending exponent scalar no 2.0

The inherited entries are elaborated in:

Options for the blending entry:

      stepwise    | Stepwise switch (discontinuous)
      max         | Maximum value switch (discontinuous)
      binomial2   | Binomial blending (smooth) n = 2
      binomial    | Binomial blending (smooth)
      exponential | Exponential blending (smooth)
      tanh        | Tanh blending (smooth)

wherein omega predictions for the viscous and inertial sublayers are blended according to the following expressions:

  • stepwise:

\[ \omega = \omega_{log} \qquad if \quad y^+ > y^+_{lam} \]

\[ \omega = \omega_{vis} \qquad if \quad y^+ <= y^+_{lam} \]

where

\( \omega \) = \(\omega\) at \(y^+\)
\( \omega_{vis} \) = \(\omega\) computed by using viscous sublayer assumptions
\( \omega_{log} \) = \(\omega\) computed by using inertial sublayer assumptions
\( y^+ \) = estimated wall-normal height of the cell centre in wall units
\( y^+_{lam} \) = estimated intersection of the viscous and inertial sublayers
  • max (PH:Eq. 27):

\[ \omega = max(\omega_{vis}, \omega_{log}) \]

  • binomial2 (ME:Eq. 15) (default):

\[ \omega = \sqrt{(\omega_{vis})^2 + (\omega_{log})^2} \]

  • binomial:

\[ \omega = ((\omega_{vis})^n + (\omega_{log})^n)^{1/n} \]

where

\( n \) = Binomial blending exponent
  • exponential (PH:Eq. 32):

\[ \omega = \omega_{vis} \exp[-\Gamma] + \omega_{log} \exp[-1/\Gamma] \]

where (PH:Eq. 31)

\( \Gamma \) = Blending expression
\( \Gamma \) = \(0.01 (y^+)^4 / (1.0 + 5.0 y^+)\)
  • tanh (KAS:Eqs. 33-34):

\[ \omega = \phi \omega_{b1} + (1 - \phi)\omega_{b2} \]

where

\( \phi \) = \(tanh((y^+/10)^4)\)
\( \omega_{b1} \) = \(\omega_{vis} + \omega_{log}\)
\( \omega_{b2} \) = \((\omega_{vis}^{1.2} + \omega_{log}^1.2)^{1/1.2}\)

G predictions for the viscous and inertial sublayers are blended in a stepwise manner, and G below \(y^+_{lam}\) (i.e. in the viscous sublayer) is presumed to be zero.

Note
  • The coefficients Cmu, kappa, and E are obtained from the specified nutWallFunction in order to ensure that each patch possesses the same set of values for these coefficients.
  • The reason why binomial2 and binomial blending methods exist at the same time is to ensure the bitwise regression with the previous versions since binomial2 and binomial with n=2 will yield slightly different output due to the miniscule differences in the implementation of the basic functions (i.e. pow, sqrt, sqr).
See also
Source files

Definition at line 291 of file omegaWallFunctionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ omegaWallFunctionFvPatchScalarField() [1/5]

omegaWallFunctionFvPatchScalarField ( const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF 
)

Construct from patch and internal field.

Definition at line 315 of file omegaWallFunctionFvPatchScalarField.C.

Referenced by omegaWallFunctionFvPatchScalarField::clone().

Here is the caller graph for this function:

◆ omegaWallFunctionFvPatchScalarField() [2/5]

omegaWallFunctionFvPatchScalarField ( const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF,
const dictionary dict 
)

Construct from patch, internal field and dictionary.

Definition at line 353 of file omegaWallFunctionFvPatchScalarField.C.

References dict, Foam::endl(), dictionary::found(), dictionary::get(), IOWarningInFunction, Foam::nl, and Foam::operator==().

Here is the call graph for this function:

◆ omegaWallFunctionFvPatchScalarField() [3/5]

omegaWallFunctionFvPatchScalarField ( const omegaWallFunctionFvPatchScalarField ptf,
const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF,
const fvPatchFieldMapper mapper 
)

Construct by mapping given omegaWallFunctionFvPatchScalarField onto a new patch

Definition at line 333 of file omegaWallFunctionFvPatchScalarField.C.

◆ omegaWallFunctionFvPatchScalarField() [4/5]

◆ omegaWallFunctionFvPatchScalarField() [5/5]

Construct as copy setting internal field reference.

Definition at line 431 of file omegaWallFunctionFvPatchScalarField.C.

◆ ~omegaWallFunctionFvPatchScalarField()

virtual ~omegaWallFunctionFvPatchScalarField ( )
virtualdefault

Destructor.

Member Function Documentation

◆ setMaster()

void setMaster ( )
protectedvirtual

Set the master patch - master is responsible for updating all wall function patches

Definition at line 55 of file omegaWallFunctionFvPatchScalarField.C.

References forAll, fvPatchField< scalar >::internalField(), omegaWallFunctionFvPatchScalarField::master(), omegaWallFunctionFvPatchScalarField::master_, omegaWallFunctionFvPatchScalarField::omega(), and omegaWallFunctionFvPatchScalarField::omegaPatch().

Here is the call graph for this function:

◆ createAveragingWeights()

void createAveragingWeights ( )
protectedvirtual

Create the averaging weights for cells which are bounded by multiple wall function faces

Definition at line 85 of file omegaWallFunctionFvPatchScalarField.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), polyMesh::changing(), Foam::dimless, forAll, mesh, IOobject::NO_READ, IOobject::NO_WRITE, fvPatchField< Type >::patchInternalField(), fvMesh::time(), Time::timeName(), and Foam::Zero.

Here is the call graph for this function:

◆ omegaPatch()

Foam::omegaWallFunctionFvPatchScalarField & omegaPatch ( const label  patchi)
protectedvirtual

Helper function to return non-const access to an omega patch.

Definition at line 145 of file omegaWallFunctionFvPatchScalarField.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryField().

Referenced by omegaWallFunctionFvPatchScalarField::setMaster().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculateTurbulenceFields()

void calculateTurbulenceFields ( const turbulenceModel turbulence,
scalarField G0,
scalarField omega0 
)
protectedvirtual

Main driver to calculate the turbulence fields.

Definition at line 162 of file omegaWallFunctionFvPatchScalarField.C.

References omegaWallFunctionFvPatchScalarField::calculate(), fvPatch::faceCells(), forAll, Foam::constant::electromagnetic::G0, and fvPatchField< Type >::patch().

Here is the call graph for this function:

◆ calculate()

void calculate ( const turbulenceModel turbulence,
const List< scalar > &  cornerWeights,
const fvPatch patch,
scalarField G,
scalarField omega 
)
protectedvirtual

Calculate the omega and G.

Reimplemented in atmOmegaWallFunctionFvPatchScalarField.

Definition at line 195 of file omegaWallFunctionFvPatchScalarField.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), nutWallFunctionFvPatchScalarField::Cmu(), Foam::exp(), forAll, Foam::constant::electromagnetic::G0, k, turbulenceModel::k(), nutWallFunctionFvPatchScalarField::kappa(), Foam::mag(), Foam::max(), turbulenceModel::nu(), nutWallFunctionFvPatchScalarField::nutw(), Foam::foamVersion::patch, Foam::pow(), Foam::pow025(), Foam::pow4(), fvPatchField< Type >::snGrad(), Foam::sqr(), Foam::sqrt(), Foam::tanh(), turbulenceModel::U(), y, turbulenceModel::y(), yPlus, and nutWallFunctionFvPatchScalarField::yPlusLam().

Referenced by omegaWallFunctionFvPatchScalarField::calculateTurbulenceFields().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ master()

virtual label& master ( )
inlineprotectedvirtual

Return non-const access to the master patch ID.

Definition at line 387 of file omegaWallFunctionFvPatchScalarField.H.

References omegaWallFunctionFvPatchScalarField::master_.

Referenced by omegaWallFunctionFvPatchScalarField::setMaster().

Here is the caller graph for this function:

◆ TypeName()

TypeName ( "omegaWallFunction"  )

Runtime type information.

◆ clone() [1/2]

virtual tmp<fvPatchScalarField> clone ( ) const
inlinevirtual

Construct and return a clone.

Reimplemented from fixedValueFvPatchField< scalar >.

Reimplemented in atmOmegaWallFunctionFvPatchScalarField.

Definition at line 434 of file omegaWallFunctionFvPatchScalarField.H.

References omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField().

Here is the call graph for this function:

◆ clone() [2/2]

virtual tmp<fvPatchScalarField> clone ( const DimensionedField< scalar, volMesh > &  iF) const
inlinevirtual

Construct and return a clone setting internal field reference.

Reimplemented from fixedValueFvPatchField< scalar >.

Reimplemented in atmOmegaWallFunctionFvPatchScalarField.

Definition at line 451 of file omegaWallFunctionFvPatchScalarField.H.

References omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField().

Here is the call graph for this function:

◆ G()

Foam::scalarField & G ( bool  init = false)

Return non-const access to the master's G field.

Definition at line 451 of file omegaWallFunctionFvPatchScalarField.C.

References init(), and Foam::foamVersion::patch.

Here is the call graph for this function:

◆ omega()

Foam::scalarField & omega ( bool  init = false)

Return non-const access to the master's omega field.

Definition at line 470 of file omegaWallFunctionFvPatchScalarField.C.

References init(), and Foam::foamVersion::patch.

Referenced by omegaWallFunctionFvPatchScalarField::setMaster().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ updateCoeffs()

void updateCoeffs ( )
virtual

Update the coefficients associated with the patch field.

Reimplemented from fvPatchField< scalar >.

Definition at line 488 of file omegaWallFunctionFvPatchScalarField.C.

References forAll, Foam::constant::universal::G, Foam::constant::electromagnetic::G0, turbulenceModel::GName(), Foam::constant::atomic::group, IOobject::groupName(), Foam::foamVersion::patch, turbulenceModel::propertiesName, and fvPatchField< Type >::updateCoeffs().

Here is the call graph for this function:

◆ updateWeightedCoeffs()

void updateWeightedCoeffs ( const scalarField weights)
virtual

Update the coefficients associated with the patch field.

Reimplemented from fvPatchField< scalar >.

Definition at line 534 of file omegaWallFunctionFvPatchScalarField.C.

References forAll, Foam::constant::universal::G, Foam::constant::electromagnetic::G0, turbulenceModel::GName(), Foam::constant::atomic::group, IOobject::groupName(), Foam::foamVersion::patch, turbulenceModel::propertiesName, and fvPatchField< Type >::updateCoeffs().

Here is the call graph for this function:

◆ manipulateMatrix() [1/2]

void manipulateMatrix ( fvMatrix< scalar > &  matrix)
virtual

Manipulate matrix.

Reimplemented from fvPatchField< scalar >.

Definition at line 591 of file omegaWallFunctionFvPatchScalarField.C.

References fvPatchField< Type >::manipulateMatrix(), Foam::foamVersion::patch, and fvMatrix< Type >::setValues().

Here is the call graph for this function:

◆ manipulateMatrix() [2/2]

void manipulateMatrix ( fvMatrix< scalar > &  matrix,
const scalarField weights 
)
virtual

Manipulate matrix with given weights.

Reimplemented from fvPatchField< scalar >.

Definition at line 607 of file omegaWallFunctionFvPatchScalarField.C.

References Foam::expressions::patchExpr::debug, Foam::endl(), fld, forAll, fvPatchField< Type >::manipulateMatrix(), Foam::foamVersion::patch, Foam::Pout, and fvMatrix< Type >::setValues().

Here is the call graph for this function:

◆ write()

void write ( Ostream os) const
virtual

Write.

Reimplemented from fixedValueFvPatchField< scalar >.

Reimplemented in atmOmegaWallFunctionFvPatchScalarField.

Definition at line 650 of file omegaWallFunctionFvPatchScalarField.C.

References os(), fixedValueFvPatchField< Type >::write(), and Ostream::writeEntry().

Here is the call graph for this function:

Member Data Documentation

◆ tolerance_

Foam::scalar tolerance_ = 1e-5
staticprotected

Tolerance used in weighted calculations.

Definition at line 327 of file omegaWallFunctionFvPatchScalarField.H.

◆ blended_

bool blended_
protected

Deprecated(2019-11) Blending switch.

Deprecated:
(2019-11) - use blending:: options

Definition at line 331 of file omegaWallFunctionFvPatchScalarField.H.

◆ initialised_

bool initialised_
protected

Initialised flag.

Definition at line 334 of file omegaWallFunctionFvPatchScalarField.H.

◆ master_

label master_
protected

◆ beta1_

scalar beta1_
protected

beta1 coefficient

Definition at line 340 of file omegaWallFunctionFvPatchScalarField.H.

◆ G_

scalarField G_
protected

Local copy of turbulence G field.

Definition at line 343 of file omegaWallFunctionFvPatchScalarField.H.

◆ omega_

scalarField omega_
protected

Local copy of turbulence omega field.

Definition at line 346 of file omegaWallFunctionFvPatchScalarField.H.

◆ cornerWeights_

List<List<scalar> > cornerWeights_
protected

List of averaging corner weights.

Definition at line 349 of file omegaWallFunctionFvPatchScalarField.H.


The documentation for this class was generated from the following files: