Go to the documentation of this file.
42 if (!isA<wallFvPatch>(
patch()))
45 <<
"Invalid wall function specification" <<
nl
46 <<
" Patch type for patch " <<
patch().
name()
47 <<
" must be wall" <<
nl
48 <<
" Current patch type is " <<
patch().type() <<
nl <<
endl
80 if (isA<epsilonWallFunctionFvPatchScalarField>(bf[
patchi]))
127 if (isA<epsilonWallFunctionFvPatchScalarField>(bf[
patchi]))
134 weights[faceCells[i]]++;
139 cornerWeights_.setSize(bf.size());
163 refCast<const epsilonWallFunctionFvPatchScalarField>(bf[
patchi]);
179 if (!cornerWeights_[
patchi].empty())
192 if (!cornerWeights_[
patchi].empty())
215 const scalar Cmu25 =
pow025(Cmu_);
216 const scalar Cmu75 =
pow(Cmu_, 0.75);
236 scalar
w = cornerWeights[facei];
238 epsilon[celli] +=
w*Cmu75*
pow(
k[celli], 1.5)/(kappa_*
y[facei]);
242 *(nutw[facei] + nuw[facei])
244 *Cmu25*
sqrt(
k[celli])
366 if (patch().index() == master_)
376 return epsilonPatch(master_).G();
385 if (patch().index() == master_)
395 return epsilonPatch(master_).epsilon(init);
417 if (patch().index() == master_)
419 createAveragingWeights();
420 calculateTurbulenceFields(turbModel,
G(
true),
epsilon(
true));
429 const_cast<FieldType&
>
431 db().lookupObject<FieldType>(turbModel.
GName())
440 G[celli] =
G0[celli];
469 if (patch().index() == master_)
471 createAveragingWeights();
472 calculateTurbulenceFields(turbModel,
G(
true),
epsilon(
true));
481 const_cast<FieldType&
>
483 db().lookupObject<FieldType>(turbModel.
GName())
493 scalar
w = weights[facei];
499 G[celli] = (1.0 -
w)*
G[celli] +
w*
G0[celli];
501 epsilonf[facei] =
epsilon[celli];
514 if (manipulatedMatrix())
519 matrix.
setValues(patch().faceCells(), patchInternalField());
531 if (manipulatedMatrix())
543 label nConstrainedCells = 0;
549 if (weights[facei] > tolerance_)
553 label celli = faceCells[facei];
555 constraintCells.
append(celli);
563 <<
": number of constrained cells = " << nConstrainedCells
564 <<
" out of " << patch().
size()
580 writeLocalEntries(os);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
bool changing() const
Is mesh changing (topology changing and/or moving)
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
static word groupName(Name name, const word &group)
Xfer< List< T > > xfer()
Transfer contents to the Xfer container as a plain List.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
virtual label & master()
Return non-const access to the master patch ID.
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))
scalar kappa_
Von Karman constant.
#define forAll(list, i)
Loop across all elements in list.
const dimensionedScalar G
Newtonian constant of gravitation.
A class for managing temporary objects.
rDeltaT dimensionedInternalField()
static scalar tolerance_
Tolerance used in weighted calculations.
static const word propertiesName
Default name of the turbulence properties dictionary.
scalar Cmu_
Cmu coefficient.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by.
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
virtual void write(Ostream &) const
Write.
dimensionedScalar pow025(const dimensionedScalar &ds)
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.
scalarField & epsilon(bool init=false)
Return non-const access to the master's epsilon field.
virtual void checkType()
Check the type of the patch.
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.
autoPtr< compressible::turbulenceModel > turbulence
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
This boundary condition provides a turbulence dissipation wall function condition for high Reynolds n...
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
virtual void write(Ostream &) const
Write.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
word GName() const
Helper function to return the name of the turbulence G field.
scalarField & G(bool init=false)
Return non-const access to the master's G field.
Abstract base class for turbulence models (RAS, LES and laminar).
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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.
Mesh data needed to do the Finite Volume discretisation.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
errorManip< error > abort(error &err)
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &epsilon0)
Main driver to calculate the turbulence fields.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual epsilonWallFunctionFvPatchScalarField & epsilonPatch(const label patchi)
Helper function to return non-const access to an epsilon patch.
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.
prefixOSstream Pout(cout, "Pout")
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
dimensionedScalar sqrt(const dimensionedScalar &ds)
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.
const dimensionedScalar e
Elementary charge.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
virtual void setMaster()
Set the master patch - master is responsible for updating all.
virtual const labelUList & faceCells() const
Return faceCells.
const Time & time() const
Return the top-level database.
Foam::fvPatchFieldMapper.
const fvPatch & patch() const
Return patch.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Generic GeometricField class.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...