Go to the documentation of this file.
76 #ifndef omegaWallFunctionFvPatchScalarField_H
77 #define omegaWallFunctionFvPatchScalarField_H
92 class omegaWallFunctionFvPatchScalarField
94 public fixedValueFvPatchField<scalar>
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by.
omegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Calculate the omega and G.
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
This boundary condition provides a wall function constraint on turbulnce specific dissipation,...
A class for managing temporary objects.
scalar kappa_
Von Karman constant.
scalarField G_
Local copy of turbulence G field.
scalar beta1_
beta1 coefficient
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
List< List< scalar > > cornerWeights_
List of averaging corner weights.
label master_
Master patch ID.
scalarField & omega(bool init=false)
Return non-const access to the master's omega field.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
scalarField omega_
Local copy of turbulence omega field.
TypeName("omegaWallFunction")
Runtime type information.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
autoPtr< compressible::turbulenceModel > turbulence
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
virtual void write(Ostream &) const
Write.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
scalar Cmu_
Cmu coefficient.
bool initialised_
Initialised flag.
Abstract base class for turbulence models (RAS, LES and laminar).
A list of keyword definitions, which are a keyword followed by any number of values (e....
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
Main driver to calculate the turbulence fields.
virtual void checkType()
Check the type of the patch.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
scalarField & G(bool init=false)
Return non-const access to the master's G field.
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
Helper function to return non-const access to an omega patch.
virtual void setMaster()
Set the master patch - master is responsible for updating all.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual label & master()
Return non-const access to the master patch ID.
Foam::fvPatchFieldMapper.
const fvPatch & patch() const
Return patch.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
static scalar tolerance_
Tolerance used in weighted calculations.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...