Go to the documentation of this file.
200 #ifndef omegaWallFunctionFvPatchScalarField_H
201 #define omegaWallFunctionFvPatchScalarField_H
216 class omegaWallFunctionFvPatchScalarField
218 public fixedValueFvPatchField<scalar>
234 static const Enum<blendingType> blendingTypeNames;
240 enum blendingType blending_;
virtual void createAveragingWeights()
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.
This boundary condition provides a wall constraint on the specific dissipation rate,...
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
A class for managing temporary objects.
virtual ~omegaWallFunctionFvPatchScalarField()=default
Destructor.
scalarField G_
Local copy of turbulence G field.
scalar beta1_
beta1 coefficient
List< List< scalar > > cornerWeights_
List of averaging corner weights.
label master_
Master patch ID.
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.
scalarField & omega(bool init=false)
Return non-const access to the master's omega field.
TypeName("omegaWallFunction")
Runtime type information.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
virtual void write(Ostream &) const
Write.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
bool initialised_
Initialised flag.
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
Helper function to return non-const access to an omega patch.
Abstract base class for turbulence models (RAS, LES and laminar).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
Main driver to calculate the turbulence fields.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
bool blended_
Deprecated(2019-11) Blending switch.
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
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...
scalarField & G(bool init=false)
Return non-const access to the master's G field.