Go to the documentation of this file.
33 #include "phaseSystem.H"
57 if (!isA<wallFvPatch>(
patch()))
60 <<
"Patch type for patch " <<
patch().name() <<
" must be wall\n"
61 <<
"Current patch type is " <<
patch().type() <<
nl
73 return 9.24*(
pow(Prat, 0.75) - 1)*(1 + 0.28*
exp(-0.007*Prat));
85 auto& ypsf = typsf.ref();
93 const scalar
f = ypt - (
log(
E_*ypt)/
kappa_ + P[facei])/Prat[facei];
94 const scalar df = 1.0 - 1.0/(ypt*
kappa_*Prat[facei]);
95 const scalar yptNew = ypt -
f/df;
103 ypsf[facei] = yptNew;
130 const label patchi =
patch().index();
133 const auto& turbModel =
164 (prevAlphat + alphaw)*hew.
snGrad()
183 auto& alphatConv = talphatConv.ref();
192 const scalar
A = qDot[facei]*rhow[facei]*
uTau[facei]*
y[facei];
193 const scalar
B = qDot[facei]*
Pr[facei]*
yPlus[facei];
200 const scalar
A = qDot[facei]*rhow[facei]*
uTau[facei]*
y[facei];
206 0.5*rhow[facei]*
uTau[facei]
212 alphatConv[facei] =
max(0.0,
alphaEff - alphaw[facei]);
247 Prt_(
dict.getOrDefault<scalar>(
"Prt", 0.85)),
248 Cmu_(
dict.getOrDefault<scalar>(
"Cmu", 0.09)),
249 kappa_(
dict.getOrDefault<scalar>(
"kappa", 0.41)),
250 E_(
dict.getOrDefault<scalar>(
"E", 9.8))
280 kappa_(awfpsf.kappa_),
295 kappa_(awfpsf.kappa_),
311 fixedValueFvPatchScalarField::updateCoeffs();
321 os.writeEntry(
"Prt",
Prt_);
322 os.writeEntry(
"Cmu",
Cmu_);
324 os.writeEntry(
"E",
E_);
326 writeEntry(
"value",
os);
335 alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
virtual tmp< Field< Type > > snGrad() const
tmp< scalarField > calcAlphat(const scalarField &prevAlphat) const
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
static const word propertiesName
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Abstract base-class for all alphatWallFunctions supporting phase-change.
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
dimensionedScalar exp(const dimensionedScalar &ds)
tmp< scalarField > Psmooth(const scalarField &Prat) const
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
dimensionedScalar pow025(const dimensionedScalar &ds)
Generic templated field type.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
virtual tmp< Field< Type > > patchInternalField() const
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
tmp< scalarField > yPlusTherm(const scalarField &P, const scalarField &Prat) const
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 write(Ostream &) const
errorManipArg< error, int > exit(error &err, const int errNo=1)
const word & name() const
void writeEntry(const word &keyword, Ostream &os) const
#define FatalErrorInFunction
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void updateCoeffs()
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)
static tmp< T > New(Args &&... args)
Foam::fvPatchFieldMapper.
Class to represent a system of phases and model interfacial transfers between them.
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.
Generic GeometricField class.
const Boundary & boundaryField() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...