Go to the documentation of this file.
50 if (!isA<wallFvPatch>(patch()))
53 <<
"Patch type for patch " << patch().name() <<
" must be wall\n"
54 <<
"Current patch type is " << patch().type() <<
nl
65 return 9.24*(
pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*
exp(-0.007*Prat));
77 for (
int i=0; i<maxIters_; i++)
79 scalar
f = ypt - (
log(E_*ypt)/kappa_ + P)/Prat;
80 scalar df = 1.0 - 1.0/(ypt*kappa_*Prat);
81 scalar yptNew = ypt -
f/df;
87 else if (
mag(yptNew - ypt) < tolerance_)
110 fixedValueFvPatchScalarField(
p, iF),
129 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
145 fixedValueFvPatchScalarField(
p, iF,
dict),
146 Prt_(
dict.lookupOrDefault<scalar>(
"Prt", 0.85)),
147 Cmu_(
dict.lookupOrDefault<scalar>(
"Cmu", 0.09)),
148 kappa_(
dict.lookupOrDefault<scalar>(
"kappa", 0.41)),
149 E_(
dict.lookupOrDefault<scalar>(
"E", 9.8))
161 fixedValueFvPatchScalarField(awfpsf),
178 fixedValueFvPatchScalarField(awfpsf, iF),
205 compressible::turbulenceModel::propertiesName,
231 turbModel.transport().he().boundaryField()[
patchi];
242 label faceCellI = patch().faceCells()[faceI];
246 scalar
yPlus =
uTau*
y[faceI]/(muw[faceI]/rhow[faceI]);
249 scalar
Pr = muw[faceI]/alphaw[faceI];
262 scalar
A = qDot[faceI]*rhow[faceI]*
uTau*
y[faceI];
263 scalar B = qDot[faceI]*
Pr*
yPlus;
269 scalar
A = qDot[faceI]*rhow[faceI]*
uTau*
y[faceI];
279 alphatw[faceI] =
max(0.0,
alphaEff - alphaw[faceI]);
284 <<
" Pr = " <<
Pr <<
nl
285 <<
" Prt = " <<
Prt_ <<
nl
286 <<
" qDot = " << qDot[faceI] <<
nl
290 <<
" alphaw = " << alphaw[faceI] <<
nl
291 <<
" alphatw = " << alphatw[faceI] <<
nl
307 writeEntry(
"value", os);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
static word groupName(Name name, const word &group)
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
scalar Prt_
Turbulent Prandtl number.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
rDeltaT dimensionedInternalField()
scalar yPlusTherm(const scalar P, const scalar Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
scalar Cmu_
Cmu coefficient.
const char *const group
Group name for atomic constants.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
simpleMatrix< scalar > A(Nc)
scalar Psmooth(const scalar Prat) const
`P' function
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
dimensionedScalar pow025(const dimensionedScalar &ds)
void write(Ostream &) const
Write.
scalar kappa_
Von Karman constant.
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.
alphatJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
A list of keyword definitions, which are a keyword followed by any number of values (e....
dimensionedScalar log(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const volScalarField & alphaEff
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void checkType()
Check the type of the patch.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< volScalarField > alphaEff() const
Return the effective turbulent thermal diffusivity for enthalpy.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
This function object evaluates and outputs turbulence y+ for turbulence models. The field is stored o...
label k
Boltzmann constant.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Graphite solid properties.
Generic GeometricField class.
virtual tmp< volScalarField > alpha() const
Return the laminar thermal diffusivity for enthalpy [kg/m/s].
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...