Go to the documentation of this file.
45 scalar alphatJayatillekeWallFunctionFvPatchScalarField::maxExp_ = 50.0;
46 scalar alphatJayatillekeWallFunctionFvPatchScalarField::tolerance_ = 0.01;
47 label alphatJayatillekeWallFunctionFvPatchScalarField::maxIters_ = 10;
51 void alphatJayatillekeWallFunctionFvPatchScalarField::checkType()
53 if (!isA<wallFvPatch>(
patch()))
56 <<
"Patch type for patch " <<
patch().name() <<
" must be wall\n"
57 <<
"Current patch type is " <<
patch().type() <<
nl
63 tmp<scalarField> alphatJayatillekeWallFunctionFvPatchScalarField::yPlus
68 const label patchi =
patch().index();
69 const tmp<volScalarField> tnut = turbModel.nut();
72 if (isA<nutWallFunctionFvPatchScalarField>(
nut.boundaryField()[patchi]))
74 const nutWallFunctionFvPatchScalarField& nutPf =
75 dynamic_cast<const nutWallFunctionFvPatchScalarField&
>
77 nut.boundaryField()[patchi]
88 y*
sqrt(turbModel.nuEff(patchi)*
mag(Uw.snGrad()))
89 /turbModel.nu(patchi);
94 scalar alphatJayatillekeWallFunctionFvPatchScalarField::Psmooth
99 return 9.24*(
pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*
exp(-0.007*Prat));
103 scalar alphatJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
111 for (
int i=0; i<maxIters_; i++)
113 scalar
f = ypt - (
log(E_*ypt)/kappa_ + P)/Prat;
114 scalar df = 1.0 - 1.0/(ypt*kappa_*Prat);
115 scalar yptNew = ypt -
f/df;
121 else if (
mag(yptNew - ypt) < tolerance_)
144 fixedValueFvPatchScalarField(
p, iF),
162 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
177 fixedValueFvPatchScalarField(
p, iF,
dict),
178 Prt_(
dict.getOrDefault<scalar>(
"Prt", 0.85)),
179 kappa_(
dict.getOrDefault<scalar>(
"kappa", 0.41)),
180 E_(
dict.getOrDefault<scalar>(
"E", 9.8))
192 fixedValueFvPatchScalarField(awfpsf),
194 kappa_(awfpsf.kappa_),
208 fixedValueFvPatchScalarField(awfpsf, iF),
210 kappa_(awfpsf.kappa_),
226 const label patchi =
patch().index();
234 compressible::turbulenceModel::propertiesName,
235 internalField().
group()
253 const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
255 turbModel.transport().he().boundaryField()[patchi];
268 scalar yPlus = yPlusp[facei];
270 scalar
uTau =
yPlus/
y[facei]*(muw[facei]/rhow[facei]);
273 scalar
Pr = muw[facei]/alphaw[facei];
276 scalar Prat =
Pr/Prt_;
279 scalar P = Psmooth(Prat);
280 scalar yPlusTherm = this->yPlusTherm(P, Prat);
284 if (yPlus < yPlusTherm)
286 scalar
A = qDot[facei]*rhow[facei]*
uTau*
y[facei];
293 scalar
A = qDot[facei]*rhow[facei]*
uTau*
y[facei];
294 scalar
B = qDot[facei]*Prt_*(1.0/kappa_*
log(E_*yPlus) + P);
295 scalar magUc =
uTau/kappa_*
log(E_*yPlusTherm) -
mag(Uw[facei]);
303 alphatw[facei] =
max(0.0,
alphaEff - alphaw[facei]);
308 <<
" Pr = " <<
Pr <<
nl
309 <<
" Prt = " << Prt_ <<
nl
310 <<
" qDot = " << qDot[facei] <<
nl
312 <<
" yPlusTherm = " << yPlusTherm <<
nl
314 <<
" alphaw = " << alphaw[facei] <<
nl
315 <<
" alphatw = " << alphatw[facei] <<
nl
327 os.writeEntry(
"Prt", Prt_);
328 os.writeEntry(
"kappa", kappa_);
329 os.writeEntry(
"E", E_);
330 writeEntry(
"value",
os);
339 alphatJayatillekeWallFunctionFvPatchScalarField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
fvPatchField< scalar > fvPatchScalarField
virtual void write(Ostream &) const
virtual void updateCoeffs()
virtual tmp< Field< Type > > snGrad() const
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr const char *const group
A class for managing temporary objects.
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Ostream & endl(Ostream &os)
dimensionedScalar exp(const dimensionedScalar &ds)
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
void write(Ostream &) const
Generic templated field type.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
alphatJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
virtual tmp< Field< Type > > patchInternalField() const
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
dimensionedScalar log(const dimensionedScalar &ds)
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
virtual void updateCoeffs()
errorManipArg< error, int > exit(error &err, const int errNo=1)
fvPatchField< vector > fvPatchVectorField
#define FatalErrorInFunction
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< volScalarField > alphaEff() const
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
Foam::fvPatchFieldMapper.
static word groupName(StringType base, const word &group)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Graphite solid properties.
virtual tmp< volScalarField > alpha() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...