Go to the documentation of this file.
48 if (!isA<wallFvPatch>(
patch()))
51 <<
"Invalid wall function specification" <<
nl
52 <<
" Patch type for patch " <<
patch().
name()
53 <<
" must be wall" <<
nl
54 <<
" Current patch type is " <<
patch().type() <<
nl <<
endl
84 if (isA<omegaWallFunctionFvPatchScalarField>(bf[
patchi]))
118 mesh.time().timeName(),
131 if (isA<omegaWallFunctionFvPatchScalarField>(bf[
patchi]))
138 label cellI = faceCells[i];
168 refCast<const omegaWallFunctionFvPatchScalarField>(bf[
patchi]);
184 if (!cornerWeights_[
patchi].empty())
197 if (!cornerWeights_[
patchi].empty())
220 const scalar Cmu25 =
pow025(Cmu_);
240 scalar
w = cornerWeights[faceI];
242 scalar omegaVis = 6.0*nuw[faceI]/(beta1_*
sqr(
y[faceI]));
244 scalar omegaLog =
sqrt(
k[cellI])/(Cmu25*kappa_*
y[faceI]);
246 omega[cellI] +=
w*
sqrt(
sqr(omegaVis) +
sqr(omegaLog));
250 *(nutw[faceI] + nuw[faceI])
252 *Cmu25*
sqrt(
k[cellI])
439 const_cast<FieldType&
>
450 G[cellI] =
G0[cellI];
451 omega[cellI] = omega0[cellI];
479 if (patch().index() == master_)
481 createAveragingWeights();
482 calculateTurbulenceFields(turbModel,
G(
true), omega(
true));
491 const_cast<FieldType&
>
493 db().lookupObject<FieldType>(turbModel.
GName())
503 scalar
w = weights[faceI];
509 G[cellI] = (1.0 -
w)*
G[cellI] +
w*
G0[cellI];
510 omega[cellI] = (1.0 -
w)*omega[cellI] +
w*omega0[cellI];
511 omegaf[faceI] = omega[cellI];
524 if (manipulatedMatrix())
529 matrix.
setValues(patch().faceCells(), patchInternalField());
541 if (manipulatedMatrix())
553 label nConstrainedCells = 0;
559 if (weights[faceI] > tolerance_)
563 label cellI = faceCells[faceI];
565 constraintCells.
append(cellI);
566 constraintomega.
append(omega[cellI]);
573 <<
": number of constrained cells = " << nConstrainedCells
574 <<
" out of " << patch().
size()
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
omegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
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.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
static word groupName(Name name, const word &group)
Xfer< List< T > > xfer()
Transfer contents to the Xfer container as a plain List.
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
This boundary condition provides a wall function constraint on turbulnce specific dissipation,...
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
#define forAll(list, i)
Loop across all elements in list.
const dimensionedScalar G
Newtonian constant of gravitation.
A class for managing temporary objects.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
const DimensionedField< Type, volMesh > & dimensionedInternalField() const
Return dimensioned internal field reference.
scalar kappa_
Von Karman constant.
rDeltaT dimensionedInternalField()
scalarField G_
Local copy of turbulence G field.
scalar beta1_
beta1 coefficient
static const word propertiesName
Default name of the turbulence properties dictionary.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const char *const group
Group name for atomic constants.
Ostream & endl(Ostream &os)
Add newline and flush stream.
List< List< scalar > > cornerWeights_
List of averaging corner weights.
label master_
Master patch ID.
const nearWallDist & y() const
Return the near wall distances.
dimensioned< scalar > mag(const dimensioned< Type > &)
bool updated() const
Return true if the boundary condition has already been updated.
scalarField & omega(bool init=false)
Return non-const access to the master's omega field.
virtual void write(Ostream &) const
Write.
dimensionedScalar pow025(const dimensionedScalar &ds)
scalarField omega_
Local copy of turbulence omega field.
const word & name() const
Return name.
virtual label size() const
Return size.
void setValues(const labelUList &cells, const UList< Type > &values)
Set solution in given cells to the specified values.
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.
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.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
bool initialised_
Initialised flag.
word GName() const
Helper function to return the name of the turbulence G field.
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....
label index() const
Return the index of this patch in the fvBoundaryMesh.
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
Main driver to calculate the turbulence fields.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Mesh data needed to do the Finite Volume discretisation.
errorManip< error > abort(error &err)
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
const double e
Elementary charge.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
prefixOSstream Pout(cout, "Pout")
const objectRegistry & db() const
Return local objectRegistry.
virtual void checkType()
Check the type of the patch.
dimensionedScalar sqrt(const dimensionedScalar &ds)
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.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
label k
Boltzmann constant.
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 const labelUList & faceCells() const
Return faceCells.
virtual label & master()
Return non-const access to the master patch ID.
Foam::fvPatchFieldMapper.
const fvPatch & patch() const
Return patch.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Generic GeometricField class.
static scalar tolerance_
Tolerance used in weighted calculations.
static scalar yPlusLam(const scalar kappa, const scalar E)
Calculate the Y+ at the edge of the laminar sublayer.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const volVectorField & U() const
Access function to velocity field.