Rationalised wall function implementation to avoid complex and inconsistent coefficients
All wall functions now operate collaboratively, obtaining the Cmu, kappa and E coefficients and yPlusLam from the nutWallFunction base class. Now these optional inputs need only be specified in the nut boundary condition with the k, epsilon, omega, v2 and f wall functions obtaining these values from there. This is much simpler to specify and avoids inconsistencies in the operation of the wall functions for the different turbulence fields. The code has also been rationalised and simplified avoiding unnecessary code and duplication.
This commit is contained in:
@ -48,9 +48,7 @@ alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField
|
||||
relax_(1.0),
|
||||
fixedDmdt_(0.0),
|
||||
L_(0.0)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField::
|
||||
|
||||
@ -24,14 +24,11 @@ License
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H"
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
#include "phaseSystem.H"
|
||||
#include "compressibleTurbulenceModel.H"
|
||||
#include "ThermalDiffusivity.H"
|
||||
#include "PhaseCompressibleTurbulenceModel.H"
|
||||
#include "wallFvPatch.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -51,18 +48,6 @@ label alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::maxIters_
|
||||
|
||||
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||
|
||||
void alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::checkType()
|
||||
{
|
||||
if (!isA<wallFvPatch>(patch()))
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Patch type for patch " << patch().name() << " must be wall\n"
|
||||
<< "Current patch type is " << patch().type() << nl
|
||||
<< exit(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
tmp<scalarField>
|
||||
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::Psmooth
|
||||
(
|
||||
@ -76,6 +61,7 @@ alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::Psmooth
|
||||
tmp<scalarField>
|
||||
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
|
||||
(
|
||||
const nutWallFunctionFvPatchScalarField& nutw,
|
||||
const scalarField& P,
|
||||
const scalarField& Prat
|
||||
) const
|
||||
@ -89,9 +75,10 @@ alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
|
||||
|
||||
for (int i=0; i<maxIters_; i++)
|
||||
{
|
||||
scalar f = ypt - (log(E_*ypt)/kappa_ + P[facei])/Prat[facei];
|
||||
scalar df = 1 - 1.0/(ypt*kappa_*Prat[facei]);
|
||||
scalar yptNew = ypt - f/df;
|
||||
const scalar f =
|
||||
ypt - (log(nutw.E()*ypt)/nutw.kappa() + P[facei])/Prat[facei];
|
||||
const scalar df = 1 - 1.0/(ypt*nutw.kappa()*Prat[facei]);
|
||||
const scalar yptNew = ypt - f/df;
|
||||
|
||||
if (yptNew < vSmall)
|
||||
{
|
||||
@ -138,7 +125,13 @@ alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::calcAlphat
|
||||
IOobject::groupName(turbulenceModel::propertiesName, phase.name())
|
||||
);
|
||||
|
||||
const scalar Cmu25 = pow025(Cmu_);
|
||||
const nutWallFunctionFvPatchScalarField& nutw =
|
||||
refCast<const nutWallFunctionFvPatchScalarField>
|
||||
(
|
||||
turbModel.nut()().boundaryField()[patchi]
|
||||
);
|
||||
|
||||
const scalar Cmu25 = pow025(nutw.Cmu());
|
||||
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
|
||||
@ -183,7 +176,7 @@ alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::calcAlphat
|
||||
// Thermal sublayer thickness
|
||||
scalarField P(this->Psmooth(Prat));
|
||||
|
||||
scalarField yPlusTherm(this->yPlusTherm(P, Prat));
|
||||
scalarField yPlusTherm(this->yPlusTherm(nutw, P, Prat));
|
||||
|
||||
tmp<scalarField> talphatConv(new scalarField(this->size()));
|
||||
scalarField& alphatConv = talphatConv.ref();
|
||||
@ -195,21 +188,31 @@ alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::calcAlphat
|
||||
scalar alphaEff = 0.0;
|
||||
if (yPlus[facei] < yPlusTherm[facei])
|
||||
{
|
||||
scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
|
||||
scalar B = qDot[facei]*Pr[facei]*yPlus[facei];
|
||||
scalar C = Pr[facei]*0.5*rhow[facei]*uTau[facei]*sqr(magUp[facei]);
|
||||
const scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
|
||||
|
||||
const scalar B = qDot[facei]*Pr[facei]*yPlus[facei];
|
||||
|
||||
const scalar C =
|
||||
Pr[facei]*0.5*rhow[facei]*uTau[facei]*sqr(magUp[facei]);
|
||||
|
||||
alphaEff = A/(B + C + vSmall);
|
||||
}
|
||||
else
|
||||
{
|
||||
scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
|
||||
scalar B =
|
||||
qDot[facei]*Prt_*(1.0/kappa_*log(E_*yPlus[facei]) + P[facei]);
|
||||
scalar magUc =
|
||||
uTau[facei]/kappa_*log(E_*yPlusTherm[facei]) - mag(Uw[facei]);
|
||||
scalar C =
|
||||
const scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
|
||||
|
||||
const scalar B =
|
||||
qDot[facei]*Prt_
|
||||
*(1.0/nutw.kappa()*log(nutw.E()*yPlus[facei]) + P[facei]);
|
||||
|
||||
const scalar magUc =
|
||||
uTau[facei]/nutw.kappa()
|
||||
*log(nutw.E()*yPlusTherm[facei]) - mag(Uw[facei]);
|
||||
|
||||
const scalar C =
|
||||
0.5*rhow[facei]*uTau[facei]
|
||||
*(Prt_*sqr(magUp[facei]) + (Pr[facei] - Prt_)*sqr(magUc));
|
||||
|
||||
alphaEff = A/(B + C + vSmall);
|
||||
}
|
||||
|
||||
@ -231,13 +234,8 @@ alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
alphatPhaseChangeWallFunctionFvPatchScalarField(p, iF),
|
||||
Prt_(0.85),
|
||||
Cmu_(0.09),
|
||||
kappa_(0.41),
|
||||
E_(9.8)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Prt_(0.85)
|
||||
{}
|
||||
|
||||
|
||||
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::
|
||||
@ -249,10 +247,7 @@ alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
alphatPhaseChangeWallFunctionFvPatchScalarField(p, iF, dict),
|
||||
Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85)),
|
||||
Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
|
||||
kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
|
||||
E_(dict.lookupOrDefault<scalar>("E", 9.8))
|
||||
Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85))
|
||||
{}
|
||||
|
||||
|
||||
@ -266,10 +261,7 @@ alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
alphatPhaseChangeWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
|
||||
Prt_(ptf.Prt_),
|
||||
Cmu_(ptf.Cmu_),
|
||||
kappa_(ptf.kappa_),
|
||||
E_(ptf.E_)
|
||||
Prt_(ptf.Prt_)
|
||||
{}
|
||||
|
||||
|
||||
@ -280,10 +272,7 @@ alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
alphatPhaseChangeWallFunctionFvPatchScalarField(awfpsf),
|
||||
Prt_(awfpsf.Prt_),
|
||||
Cmu_(awfpsf.Cmu_),
|
||||
kappa_(awfpsf.kappa_),
|
||||
E_(awfpsf.E_)
|
||||
Prt_(awfpsf.Prt_)
|
||||
{}
|
||||
|
||||
|
||||
@ -295,10 +284,7 @@ alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
alphatPhaseChangeWallFunctionFvPatchScalarField(awfpsf, iF),
|
||||
Prt_(awfpsf.Prt_),
|
||||
Cmu_(awfpsf.Cmu_),
|
||||
kappa_(awfpsf.kappa_),
|
||||
E_(awfpsf.E_)
|
||||
Prt_(awfpsf.Prt_)
|
||||
{}
|
||||
|
||||
|
||||
@ -324,9 +310,6 @@ void alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::write
|
||||
{
|
||||
fvPatchField<scalar>::write(os);
|
||||
writeEntry(os, "Prt", Prt_);
|
||||
writeEntry(os, "Cmu", Cmu_);
|
||||
writeEntry(os, "kappa", kappa_);
|
||||
writeEntry(os, "E", E_);
|
||||
writeEntry(os, "dmdt", dmdt_);
|
||||
writeEntry(os, "value", *this);
|
||||
}
|
||||
|
||||
@ -63,6 +63,7 @@ SourceFiles
|
||||
#define compressible_alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField_H
|
||||
|
||||
#include "alphatPhaseChangeWallFunctionFvPatchScalarField.H"
|
||||
#include "nutWallFunctionFvPatchScalarField.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -87,15 +88,6 @@ protected:
|
||||
//- Turbulent Prandtl number
|
||||
scalar Prt_;
|
||||
|
||||
//- Cmu coefficient
|
||||
scalar Cmu_;
|
||||
|
||||
//- Von Karman constant
|
||||
scalar kappa_;
|
||||
|
||||
//- E coefficient
|
||||
scalar E_;
|
||||
|
||||
// Solution parameters
|
||||
|
||||
static scalar maxExp_;
|
||||
@ -105,15 +97,13 @@ protected:
|
||||
|
||||
// Protected Member Functions
|
||||
|
||||
//- Check the type of the patch
|
||||
void checkType();
|
||||
|
||||
//- 'P' function
|
||||
tmp<scalarField> Psmooth(const scalarField& Prat) const;
|
||||
|
||||
//- Calculate y+ at the edge of the thermal laminar sublayer
|
||||
tmp<scalarField> yPlusTherm
|
||||
(
|
||||
const nutWallFunctionFvPatchScalarField& nutw,
|
||||
const scalarField& P,
|
||||
const scalarField& Prat
|
||||
) const;
|
||||
|
||||
@ -24,17 +24,12 @@ License
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "alphatWallBoilingWallFunctionFvPatchScalarField.H"
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
#include "phaseSystem.H"
|
||||
#include "compressibleTurbulenceModel.H"
|
||||
#include "ThermalDiffusivity.H"
|
||||
#include "PhaseCompressibleTurbulenceModel.H"
|
||||
#include "saturationModel.H"
|
||||
#include "wallFvPatch.H"
|
||||
#include "uniformDimensionedFields.H"
|
||||
#include "mathematicalConstants.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
using namespace Foam::constant::mathematical;
|
||||
|
||||
@ -428,9 +423,13 @@ void alphatWallBoilingWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
)
|
||||
);
|
||||
|
||||
const tmp<scalarField> tnutw = turbModel.nut(patchi);
|
||||
const nutWallFunctionFvPatchScalarField& nutw =
|
||||
refCast<const nutWallFunctionFvPatchScalarField>
|
||||
(
|
||||
turbModel.nut()().boundaryField()[patchi]
|
||||
);
|
||||
|
||||
const scalar Cmu25(pow025(Cmu_));
|
||||
const scalar Cmu25(pow025(nutw.Cmu()));
|
||||
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
|
||||
@ -474,7 +473,7 @@ void alphatWallBoilingWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
// Thermal sublayer thickness
|
||||
const scalarField P(this->Psmooth(Prat));
|
||||
|
||||
const scalarField yPlusTherm(this->yPlusTherm(P, Prat));
|
||||
const scalarField yPlusTherm(this->yPlusTherm(nutw, P, Prat));
|
||||
|
||||
tmp<volScalarField> tCp = liquid.thermo().Cp();
|
||||
const volScalarField& Cp = tCp();
|
||||
@ -517,11 +516,24 @@ void alphatWallBoilingWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
{
|
||||
// Liquid temperature at y+=250 is estimated from logarithmic
|
||||
// thermal wall function (Koncar, Krepper & Egorov, 2005)
|
||||
const scalarField Tplus_y250(Prt_*(log(E_*250)/kappa_ + P));
|
||||
const scalarField Tplus_y250
|
||||
(
|
||||
Prt_*(log(nutw.E()*250)/nutw.kappa() + P)
|
||||
);
|
||||
|
||||
const scalarField Tplus(Prt_*(log(E_*yPlus)/kappa_ + P));
|
||||
scalarField Tl(Tw - (Tplus_y250/Tplus)*(Tw - Tc));
|
||||
Tl = max(Tc - 40, Tl);
|
||||
const scalarField Tplus
|
||||
(
|
||||
Prt_*(log(nutw.E()*yPlus)/nutw.kappa() + P)
|
||||
);
|
||||
|
||||
const scalarField Tl
|
||||
(
|
||||
max
|
||||
(
|
||||
Tc - 40,
|
||||
Tw - (Tplus_y250/Tplus)*(Tw - Tc)
|
||||
)
|
||||
);
|
||||
|
||||
// Nucleation site density:
|
||||
const scalarField N
|
||||
|
||||
@ -25,10 +25,7 @@ License
|
||||
|
||||
#include "alphatJayatillekeWallFunctionFvPatchScalarField.H"
|
||||
#include "turbulentFluidThermoModel.H"
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "volFields.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
#include "wallFvPatch.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -45,18 +42,6 @@ label alphatJayatillekeWallFunctionFvPatchScalarField::maxIters_ = 10;
|
||||
|
||||
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||
|
||||
void alphatJayatillekeWallFunctionFvPatchScalarField::checkType()
|
||||
{
|
||||
if (!isA<wallFvPatch>(patch()))
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Patch type for patch " << patch().name() << " must be wall\n"
|
||||
<< "Current patch type is " << patch().type() << nl
|
||||
<< exit(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
scalar alphatJayatillekeWallFunctionFvPatchScalarField::Psmooth
|
||||
(
|
||||
const scalar Prat
|
||||
@ -68,6 +53,7 @@ scalar alphatJayatillekeWallFunctionFvPatchScalarField::Psmooth
|
||||
|
||||
scalar alphatJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
|
||||
(
|
||||
const nutWallFunctionFvPatchScalarField& nutw,
|
||||
const scalar P,
|
||||
const scalar Prat
|
||||
) const
|
||||
@ -76,8 +62,8 @@ scalar alphatJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
|
||||
|
||||
for (int i=0; i<maxIters_; i++)
|
||||
{
|
||||
scalar f = ypt - (log(E_*ypt)/kappa_ + P)/Prat;
|
||||
scalar df = 1.0 - 1.0/(ypt*kappa_*Prat);
|
||||
scalar f = ypt - (log(nutw.E()*ypt)/nutw.kappa() + P)/Prat;
|
||||
scalar df = 1.0 - 1.0/(ypt*nutw.kappa()*Prat);
|
||||
scalar yptNew = ypt - f/df;
|
||||
|
||||
if (yptNew < vSmall)
|
||||
@ -108,13 +94,8 @@ alphatJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(p, iF),
|
||||
Prt_(0.85),
|
||||
Cmu_(0.09),
|
||||
kappa_(0.41),
|
||||
E_(9.8)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Prt_(0.85)
|
||||
{}
|
||||
|
||||
|
||||
alphatJayatillekeWallFunctionFvPatchScalarField::
|
||||
@ -127,10 +108,7 @@ alphatJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(ptf, p, iF, mapper),
|
||||
Prt_(ptf.Prt_),
|
||||
Cmu_(ptf.Cmu_),
|
||||
kappa_(ptf.kappa_),
|
||||
E_(ptf.E_)
|
||||
Prt_(ptf.Prt_)
|
||||
{}
|
||||
|
||||
|
||||
@ -143,13 +121,8 @@ alphatJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(p, iF, dict),
|
||||
Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85)),
|
||||
Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
|
||||
kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
|
||||
E_(dict.lookupOrDefault<scalar>("E", 9.8))
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85))
|
||||
{}
|
||||
|
||||
|
||||
alphatJayatillekeWallFunctionFvPatchScalarField::
|
||||
@ -159,13 +132,8 @@ alphatJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(awfpsf),
|
||||
Prt_(awfpsf.Prt_),
|
||||
Cmu_(awfpsf.Cmu_),
|
||||
kappa_(awfpsf.kappa_),
|
||||
E_(awfpsf.E_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Prt_(awfpsf.Prt_)
|
||||
{}
|
||||
|
||||
|
||||
alphatJayatillekeWallFunctionFvPatchScalarField::
|
||||
@ -176,13 +144,8 @@ alphatJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(awfpsf, iF),
|
||||
Prt_(awfpsf.Prt_),
|
||||
Cmu_(awfpsf.Cmu_),
|
||||
kappa_(awfpsf.kappa_),
|
||||
E_(awfpsf.E_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Prt_(awfpsf.Prt_)
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
@ -207,7 +170,13 @@ void alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
)
|
||||
);
|
||||
|
||||
const scalar Cmu25 = pow025(Cmu_);
|
||||
const nutWallFunctionFvPatchScalarField& nutw =
|
||||
refCast<const nutWallFunctionFvPatchScalarField>
|
||||
(
|
||||
turbModel.nut()().boundaryField()[patchi]
|
||||
);
|
||||
|
||||
const scalar Cmu25 = pow025(nutw.Cmu());
|
||||
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
|
||||
@ -253,25 +222,34 @@ void alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
|
||||
// Thermal sublayer thickness
|
||||
scalar P = Psmooth(Prat);
|
||||
scalar yPlusTherm = this->yPlusTherm(P, Prat);
|
||||
scalar yPlusTherm = this->yPlusTherm(nutw, P, Prat);
|
||||
|
||||
// Evaluate new effective thermal diffusivity
|
||||
scalar alphaEff = 0.0;
|
||||
if (yPlus < yPlusTherm)
|
||||
{
|
||||
scalar A = qDot[facei]*rhow[facei]*uTau*y[facei];
|
||||
scalar B = qDot[facei]*Pr*yPlus;
|
||||
scalar C = Pr*0.5*rhow[facei]*uTau*sqr(magUp[facei]);
|
||||
const scalar A = qDot[facei]*rhow[facei]*uTau*y[facei];
|
||||
|
||||
const scalar B = qDot[facei]*Pr*yPlus;
|
||||
|
||||
const scalar C = Pr*0.5*rhow[facei]*uTau*sqr(magUp[facei]);
|
||||
|
||||
alphaEff = A/(B + C + vSmall);
|
||||
}
|
||||
else
|
||||
{
|
||||
scalar A = qDot[facei]*rhow[facei]*uTau*y[facei];
|
||||
scalar B = qDot[facei]*Prt_*(1.0/kappa_*log(E_*yPlus) + P);
|
||||
scalar magUc = uTau/kappa_*log(E_*yPlusTherm) - mag(Uw[facei]);
|
||||
scalar C =
|
||||
const scalar A = qDot[facei]*rhow[facei]*uTau*y[facei];
|
||||
|
||||
const scalar B =
|
||||
qDot[facei]*Prt_*(1.0/nutw.kappa()*log(nutw.E()*yPlus) + P);
|
||||
|
||||
const scalar magUc =
|
||||
uTau/nutw.kappa()*log(nutw.E()*yPlusTherm) - mag(Uw[facei]);
|
||||
|
||||
const scalar C =
|
||||
0.5*rhow[facei]*uTau
|
||||
*(Prt_*sqr(magUp[facei]) + (Pr - Prt_)*sqr(magUc));
|
||||
|
||||
alphaEff = A/(B + C + vSmall);
|
||||
}
|
||||
|
||||
@ -301,9 +279,6 @@ void alphatJayatillekeWallFunctionFvPatchScalarField::write(Ostream& os) const
|
||||
{
|
||||
fvPatchField<scalar>::write(os);
|
||||
writeEntry(os, "Prt", Prt_);
|
||||
writeEntry(os, "Cmu", Cmu_);
|
||||
writeEntry(os, "kappa", kappa_);
|
||||
writeEntry(os, "E", E_);
|
||||
writeEntry(os, "value", *this);
|
||||
}
|
||||
|
||||
|
||||
@ -61,6 +61,7 @@ SourceFiles
|
||||
#define alphatJayatillekeWallFunctionFvPatchScalarField_H
|
||||
|
||||
#include "fixedValueFvPatchFields.H"
|
||||
#include "nutWallFunctionFvPatchScalarField.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -82,15 +83,6 @@ class alphatJayatillekeWallFunctionFvPatchScalarField
|
||||
//- Turbulent Prandtl number
|
||||
scalar Prt_;
|
||||
|
||||
//- Cmu coefficient
|
||||
scalar Cmu_;
|
||||
|
||||
//- Von Karman constant
|
||||
scalar kappa_;
|
||||
|
||||
//- E coefficient
|
||||
scalar E_;
|
||||
|
||||
|
||||
// Solution parameters
|
||||
|
||||
@ -101,15 +93,13 @@ class alphatJayatillekeWallFunctionFvPatchScalarField
|
||||
|
||||
// Private Member Functions
|
||||
|
||||
//- Check the type of the patch
|
||||
void checkType();
|
||||
|
||||
//- `P' function
|
||||
scalar Psmooth(const scalar Prat) const;
|
||||
|
||||
//- Calculate y+ at the edge of the thermal laminar sublayer
|
||||
scalar yPlusTherm
|
||||
(
|
||||
const nutWallFunctionFvPatchScalarField& nutw,
|
||||
const scalar P,
|
||||
const scalar Prat
|
||||
) const;
|
||||
|
||||
@ -25,9 +25,6 @@ License
|
||||
|
||||
#include "alphatJayatillekeWallFunctionFvPatchScalarField.H"
|
||||
#include "turbulenceModel.H"
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "volFields.H"
|
||||
#include "wallFvPatch.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
@ -44,20 +41,6 @@ label alphatJayatillekeWallFunctionFvPatchScalarField::maxIters_ = 10;
|
||||
|
||||
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
|
||||
|
||||
void alphatJayatillekeWallFunctionFvPatchScalarField::checkType()
|
||||
{
|
||||
if (!isA<wallFvPatch>(patch()))
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Invalid wall function specification" << nl
|
||||
<< " Patch type for patch " << patch().name()
|
||||
<< " must be wall" << nl
|
||||
<< " Current patch type is " << patch().type() << nl << endl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
scalar alphatJayatillekeWallFunctionFvPatchScalarField::Psmooth
|
||||
(
|
||||
const scalar Prat
|
||||
@ -69,6 +52,7 @@ scalar alphatJayatillekeWallFunctionFvPatchScalarField::Psmooth
|
||||
|
||||
scalar alphatJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
|
||||
(
|
||||
const nutWallFunctionFvPatchScalarField& nutw,
|
||||
const scalar P,
|
||||
const scalar Prat
|
||||
) const
|
||||
@ -77,8 +61,8 @@ scalar alphatJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
|
||||
|
||||
for (int i=0; i<maxIters_; i++)
|
||||
{
|
||||
scalar f = ypt - (log(E_*ypt)/kappa_ + P)/Prat;
|
||||
scalar df = 1.0 - 1.0/(ypt*kappa_*Prat);
|
||||
scalar f = ypt - (log(nutw.E()*ypt)/nutw.kappa() + P)/Prat;
|
||||
scalar df = 1.0 - 1.0/(ypt*nutw.kappa()*Prat);
|
||||
scalar yptNew = ypt - f/df;
|
||||
|
||||
if (yptNew < vSmall)
|
||||
@ -109,13 +93,8 @@ alphatJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(p, iF),
|
||||
Prt_(0.85),
|
||||
Cmu_(0.09),
|
||||
kappa_(0.41),
|
||||
E_(9.8)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Prt_(0.85)
|
||||
{}
|
||||
|
||||
|
||||
alphatJayatillekeWallFunctionFvPatchScalarField::
|
||||
@ -128,13 +107,8 @@ alphatJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(ptf, p, iF, mapper),
|
||||
Prt_(ptf.Prt_),
|
||||
Cmu_(ptf.Cmu_),
|
||||
kappa_(ptf.kappa_),
|
||||
E_(ptf.E_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Prt_(ptf.Prt_)
|
||||
{}
|
||||
|
||||
|
||||
alphatJayatillekeWallFunctionFvPatchScalarField::
|
||||
@ -146,13 +120,8 @@ alphatJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(p, iF, dict),
|
||||
Prt_(readScalar(dict.lookup("Prt"))), // force read to avoid ambiguity
|
||||
Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
|
||||
kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
|
||||
E_(dict.lookupOrDefault<scalar>("E", 9.8))
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Prt_(readScalar(dict.lookup("Prt"))) // force read to avoid ambiguity
|
||||
{}
|
||||
|
||||
|
||||
alphatJayatillekeWallFunctionFvPatchScalarField::
|
||||
@ -162,13 +131,8 @@ alphatJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(wfpsf),
|
||||
Prt_(wfpsf.Prt_),
|
||||
Cmu_(wfpsf.Cmu_),
|
||||
kappa_(wfpsf.kappa_),
|
||||
E_(wfpsf.E_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Prt_(wfpsf.Prt_)
|
||||
{}
|
||||
|
||||
|
||||
alphatJayatillekeWallFunctionFvPatchScalarField::
|
||||
@ -179,13 +143,8 @@ alphatJayatillekeWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(wfpsf, iF),
|
||||
Prt_(wfpsf.Prt_),
|
||||
Cmu_(wfpsf.Cmu_),
|
||||
kappa_(wfpsf.kappa_),
|
||||
E_(wfpsf.E_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Prt_(wfpsf.Prt_)
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
@ -210,7 +169,13 @@ void alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
)
|
||||
);
|
||||
|
||||
const scalar Cmu25 = pow(Cmu_, 0.25);
|
||||
const nutWallFunctionFvPatchScalarField& nutw =
|
||||
refCast<const nutWallFunctionFvPatchScalarField>
|
||||
(
|
||||
turbModel.nut()().boundaryField()[patchi]
|
||||
);
|
||||
|
||||
const scalar Cmu25 = pow(nutw.Cmu(), 0.25);
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
const tmp<volScalarField> tnu = turbModel.nu();
|
||||
const volScalarField& nu = tnu();
|
||||
@ -236,23 +201,24 @@ void alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
scalarField& alphatw = *this;
|
||||
forAll(alphatw, facei)
|
||||
{
|
||||
label celli = patch().faceCells()[facei];
|
||||
const label celli = patch().faceCells()[facei];
|
||||
|
||||
// y+
|
||||
scalar yPlus = Cmu25*sqrt(k[celli])*y[facei]/nuw[facei];
|
||||
const scalar yPlus = Cmu25*sqrt(k[celli])*y[facei]/nuw[facei];
|
||||
|
||||
// Molecular-to-turbulent Prandtl number ratio
|
||||
scalar Prat = Pr/Prt_;
|
||||
const scalar Prat = Pr/Prt_;
|
||||
|
||||
// Thermal sublayer thickness
|
||||
scalar P = Psmooth(Prat);
|
||||
scalar yPlusTherm = this->yPlusTherm(P, Prat);
|
||||
const scalar P = Psmooth(Prat);
|
||||
const scalar yPlusTherm = this->yPlusTherm(nutw, P, Prat);
|
||||
|
||||
// Update turbulent thermal conductivity
|
||||
if (yPlus > yPlusTherm)
|
||||
{
|
||||
scalar nu = nuw[facei];
|
||||
scalar kt = nu*(yPlus/(Prt_*(log(E_*yPlus)/kappa_ + P)) - 1/Pr);
|
||||
const scalar nu = nuw[facei];
|
||||
const scalar kt =
|
||||
nu*(yPlus/(Prt_*(log(nutw.E()*yPlus)/nutw.kappa() + P)) - 1/Pr);
|
||||
|
||||
alphatw[facei] = max(0.0, kt);
|
||||
}
|
||||
else
|
||||
@ -269,9 +235,6 @@ void alphatJayatillekeWallFunctionFvPatchScalarField::write(Ostream& os) const
|
||||
{
|
||||
fvPatchField<scalar>::write(os);
|
||||
writeEntry(os, "Prt", Prt_);
|
||||
writeEntry(os, "Cmu", Cmu_);
|
||||
writeEntry(os, "kappa", kappa_);
|
||||
writeEntry(os, "E", E_);
|
||||
writeEntry(os, "value", *this);
|
||||
}
|
||||
|
||||
|
||||
@ -60,6 +60,7 @@ SourceFiles
|
||||
#define alphatJayatillekeWallFunctionFvPatchScalarField_H
|
||||
|
||||
#include "fixedValueFvPatchFields.H"
|
||||
#include "nutWallFunctionFvPatchScalarField.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -101,15 +102,13 @@ protected:
|
||||
|
||||
// Protected member functions
|
||||
|
||||
//- Check the type of the patch
|
||||
virtual void checkType();
|
||||
|
||||
//- `P' function
|
||||
scalar Psmooth(const scalar Prat) const;
|
||||
|
||||
//- Calculate y+ at the edge of the thermal laminar sublayer
|
||||
scalar yPlusTherm
|
||||
(
|
||||
const nutWallFunctionFvPatchScalarField& nutw,
|
||||
const scalar P,
|
||||
const scalar Prat
|
||||
) const;
|
||||
|
||||
@ -26,10 +26,7 @@ License
|
||||
#include "epsilonWallFunctionFvPatchScalarField.H"
|
||||
#include "nutWallFunctionFvPatchScalarField.H"
|
||||
#include "turbulenceModel.H"
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "fvMatrix.H"
|
||||
#include "volFields.H"
|
||||
#include "wallFvPatch.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
@ -38,31 +35,6 @@ Foam::scalar Foam::epsilonWallFunctionFvPatchScalarField::tolerance_ = 1e-5;
|
||||
|
||||
// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
|
||||
|
||||
void Foam::epsilonWallFunctionFvPatchScalarField::checkType()
|
||||
{
|
||||
if (!isA<wallFvPatch>(patch()))
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Invalid wall function specification" << nl
|
||||
<< " Patch type for patch " << patch().name()
|
||||
<< " must be wall" << nl
|
||||
<< " Current patch type is " << patch().type() << nl << endl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
void Foam::epsilonWallFunctionFvPatchScalarField::writeLocalEntries
|
||||
(
|
||||
Ostream& os
|
||||
) const
|
||||
{
|
||||
writeEntry(os, "Cmu", Cmu_);
|
||||
writeEntry(os, "kappa", kappa_);
|
||||
writeEntry(os, "E", E_);
|
||||
}
|
||||
|
||||
|
||||
void Foam::epsilonWallFunctionFvPatchScalarField::setMaster()
|
||||
{
|
||||
if (master_ != -1)
|
||||
@ -211,24 +183,27 @@ void Foam::epsilonWallFunctionFvPatchScalarField::calculate
|
||||
{
|
||||
const label patchi = patch.index();
|
||||
|
||||
const nutWallFunctionFvPatchScalarField& nutw =
|
||||
refCast<const nutWallFunctionFvPatchScalarField>
|
||||
(
|
||||
turbModel.nut()().boundaryField()[patchi]
|
||||
);
|
||||
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
|
||||
const scalar Cmu25 = pow025(Cmu_);
|
||||
const scalar Cmu75 = pow(Cmu_, 0.75);
|
||||
|
||||
const tmp<volScalarField> tk = turbModel.k();
|
||||
const volScalarField& k = tk();
|
||||
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
|
||||
const tmp<scalarField> tnutw = turbModel.nut(patchi);
|
||||
const scalarField& nutw = tnutw();
|
||||
const tmp<volScalarField> tk = turbModel.k();
|
||||
const volScalarField& k = tk();
|
||||
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
|
||||
const scalarField magGradUw(mag(Uw.snGrad()));
|
||||
|
||||
const scalar Cmu25 = pow025(nutw.Cmu());
|
||||
const scalar Cmu75 = pow(nutw.Cmu(), 0.75);
|
||||
|
||||
// Set epsilon and G
|
||||
forAll(nutw, facei)
|
||||
{
|
||||
@ -238,16 +213,17 @@ void Foam::epsilonWallFunctionFvPatchScalarField::calculate
|
||||
|
||||
const scalar w = cornerWeights[facei];
|
||||
|
||||
if (yPlus > yPlusLam_)
|
||||
if (yPlus > nutw.yPlusLam())
|
||||
{
|
||||
epsilon0[celli] += w*Cmu75*pow(k[celli], 1.5)/(kappa_*y[facei]);
|
||||
epsilon0[celli] +=
|
||||
w*Cmu75*pow(k[celli], 1.5)/(nutw.kappa()*y[facei]);
|
||||
|
||||
G0[celli] +=
|
||||
w
|
||||
*(nutw[facei] + nuw[facei])
|
||||
*magGradUw[facei]
|
||||
*Cmu25*sqrt(k[celli])
|
||||
/(kappa_*y[facei]);
|
||||
/(nutw.kappa()*y[facei]);
|
||||
}
|
||||
else
|
||||
{
|
||||
@ -267,18 +243,12 @@ epsilonWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(p, iF),
|
||||
Cmu_(0.09),
|
||||
kappa_(0.41),
|
||||
E_(9.8),
|
||||
yPlusLam_(nutWallFunctionFvPatchScalarField::yPlusLam(kappa_, E_)),
|
||||
G_(),
|
||||
epsilon_(),
|
||||
initialised_(false),
|
||||
master_(-1),
|
||||
cornerWeights_()
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
Foam::epsilonWallFunctionFvPatchScalarField::
|
||||
@ -290,18 +260,12 @@ epsilonWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(p, iF, dict),
|
||||
Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
|
||||
kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
|
||||
E_(dict.lookupOrDefault<scalar>("E", 9.8)),
|
||||
yPlusLam_(nutWallFunctionFvPatchScalarField::yPlusLam(kappa_, E_)),
|
||||
G_(),
|
||||
epsilon_(),
|
||||
initialised_(false),
|
||||
master_(-1),
|
||||
cornerWeights_()
|
||||
{
|
||||
checkType();
|
||||
|
||||
// Apply zero-gradient condition on start-up
|
||||
this->operator==(patchInternalField());
|
||||
}
|
||||
@ -317,18 +281,12 @@ epsilonWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
|
||||
Cmu_(ptf.Cmu_),
|
||||
kappa_(ptf.kappa_),
|
||||
E_(ptf.E_),
|
||||
yPlusLam_(ptf.yPlusLam_),
|
||||
G_(),
|
||||
epsilon_(),
|
||||
initialised_(false),
|
||||
master_(-1),
|
||||
cornerWeights_()
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
Foam::epsilonWallFunctionFvPatchScalarField::
|
||||
@ -338,18 +296,12 @@ epsilonWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(ewfpsf),
|
||||
Cmu_(ewfpsf.Cmu_),
|
||||
kappa_(ewfpsf.kappa_),
|
||||
E_(ewfpsf.E_),
|
||||
yPlusLam_(ewfpsf.yPlusLam_),
|
||||
G_(),
|
||||
epsilon_(),
|
||||
initialised_(false),
|
||||
master_(-1),
|
||||
cornerWeights_()
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
Foam::epsilonWallFunctionFvPatchScalarField::
|
||||
@ -360,18 +312,12 @@ epsilonWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(ewfpsf, iF),
|
||||
Cmu_(ewfpsf.Cmu_),
|
||||
kappa_(ewfpsf.kappa_),
|
||||
E_(ewfpsf.E_),
|
||||
yPlusLam_(ewfpsf.yPlusLam_),
|
||||
G_(),
|
||||
epsilon_(),
|
||||
initialised_(false),
|
||||
master_(-1),
|
||||
cornerWeights_()
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
@ -590,13 +536,6 @@ void Foam::epsilonWallFunctionFvPatchScalarField::manipulateMatrix
|
||||
}
|
||||
|
||||
|
||||
void Foam::epsilonWallFunctionFvPatchScalarField::write(Ostream& os) const
|
||||
{
|
||||
writeLocalEntries(os);
|
||||
fixedValueFvPatchField<scalar>::write(os);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
|
||||
@ -40,16 +40,10 @@ Description
|
||||
\endvartable
|
||||
|
||||
The model switches between laminar and turbulent functions based on the
|
||||
laminar-to-turbulent y+ value derived from kappa and E.
|
||||
laminar-to-turbulent y+ value derived from the kappa and E specified in the
|
||||
corresponding nutWallFunction.
|
||||
|
||||
Usage
|
||||
\table
|
||||
Property | Description | Required | Default value
|
||||
Cmu | model coefficient | no | 0.09
|
||||
kappa | Von Karman constant | no | 0.41
|
||||
E | model coefficient | no | 9.8
|
||||
\endtable
|
||||
|
||||
Example of the boundary condition specification:
|
||||
\verbatim
|
||||
<patchName>
|
||||
@ -94,18 +88,6 @@ protected:
|
||||
//- Tolerance used in weighted calculations
|
||||
static scalar tolerance_;
|
||||
|
||||
//- Cmu coefficient
|
||||
scalar Cmu_;
|
||||
|
||||
//- Von Karman constant
|
||||
scalar kappa_;
|
||||
|
||||
//- E coefficient
|
||||
scalar E_;
|
||||
|
||||
//- y+ at the edge of the laminar sublayer
|
||||
scalar yPlusLam_;
|
||||
|
||||
//- Local copy of turbulence G field
|
||||
scalarField G_;
|
||||
|
||||
@ -124,12 +106,6 @@ protected:
|
||||
|
||||
// Protected Member Functions
|
||||
|
||||
//- Check the type of the patch
|
||||
virtual void checkType();
|
||||
|
||||
//- Write local wall function variables
|
||||
virtual void writeLocalEntries(Ostream&) const;
|
||||
|
||||
//- Set the master patch - master is responsible for updating all
|
||||
// wall function patches
|
||||
virtual void setMaster();
|
||||
@ -270,12 +246,6 @@ public:
|
||||
fvMatrix<scalar>& matrix,
|
||||
const scalarField& weights
|
||||
);
|
||||
|
||||
|
||||
// I-O
|
||||
|
||||
//- Write
|
||||
virtual void write(Ostream&) const;
|
||||
};
|
||||
|
||||
|
||||
|
||||
@ -24,9 +24,8 @@ License
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "fWallFunctionFvPatchScalarField.H"
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "volFields.H"
|
||||
#include "wallFvPatch.H"
|
||||
#include "nutWallFunctionFvPatchScalarField.H"
|
||||
#include "turbulenceModel.H"
|
||||
#include "v2f.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
@ -37,47 +36,6 @@ namespace Foam
|
||||
namespace RASModels
|
||||
{
|
||||
|
||||
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
|
||||
|
||||
void fWallFunctionFvPatchScalarField::checkType()
|
||||
{
|
||||
if (!isA<wallFvPatch>(patch()))
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Invalid wall function specification" << nl
|
||||
<< " Patch type for patch " << patch().name()
|
||||
<< " must be wall" << nl
|
||||
<< " Current patch type is " << patch().type() << nl << endl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
void fWallFunctionFvPatchScalarField::writeLocalEntries(Ostream& os) const
|
||||
{
|
||||
writeEntry(os, "Cmu", Cmu_);
|
||||
writeEntry(os, "kappa", kappa_);
|
||||
writeEntry(os, "E", E_);
|
||||
}
|
||||
|
||||
|
||||
scalar fWallFunctionFvPatchScalarField::yPlusLam
|
||||
(
|
||||
const scalar kappa,
|
||||
const scalar E
|
||||
)
|
||||
{
|
||||
scalar ypl = 11.0;
|
||||
|
||||
for (int i=0; i<10; i++)
|
||||
{
|
||||
ypl = log(max(E*ypl, 1))/kappa;
|
||||
}
|
||||
|
||||
return ypl;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField
|
||||
@ -86,14 +44,8 @@ fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField
|
||||
const DimensionedField<scalar, volMesh>& iF
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(p, iF),
|
||||
Cmu_(0.09),
|
||||
kappa_(0.41),
|
||||
E_(9.8),
|
||||
yPlusLam_(yPlusLam(kappa_, E_))
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
fixedValueFvPatchField<scalar>(p, iF)
|
||||
{}
|
||||
|
||||
|
||||
fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField
|
||||
@ -104,14 +56,8 @@ fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField
|
||||
const fvPatchFieldMapper& mapper
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
|
||||
Cmu_(ptf.Cmu_),
|
||||
kappa_(ptf.kappa_),
|
||||
E_(ptf.E_),
|
||||
yPlusLam_(ptf.yPlusLam_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
fixedValueFvPatchField<scalar>(ptf, p, iF, mapper)
|
||||
{}
|
||||
|
||||
|
||||
fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField
|
||||
@ -121,14 +67,8 @@ fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField
|
||||
const dictionary& dict
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(p, iF, dict),
|
||||
Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
|
||||
kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
|
||||
E_(dict.lookupOrDefault<scalar>("E", 9.8)),
|
||||
yPlusLam_(yPlusLam(kappa_, E_))
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
fixedValueFvPatchField<scalar>(p, iF, dict)
|
||||
{}
|
||||
|
||||
|
||||
fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField
|
||||
@ -136,14 +76,8 @@ fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField
|
||||
const fWallFunctionFvPatchScalarField& v2wfpsf
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(v2wfpsf),
|
||||
Cmu_(v2wfpsf.Cmu_),
|
||||
kappa_(v2wfpsf.kappa_),
|
||||
E_(v2wfpsf.E_),
|
||||
yPlusLam_(v2wfpsf.yPlusLam_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
fixedValueFvPatchField<scalar>(v2wfpsf)
|
||||
{}
|
||||
|
||||
|
||||
fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField
|
||||
@ -152,14 +86,8 @@ fWallFunctionFvPatchScalarField::fWallFunctionFvPatchScalarField
|
||||
const DimensionedField<scalar, volMesh>& iF
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(v2wfpsf, iF),
|
||||
Cmu_(v2wfpsf.Cmu_),
|
||||
kappa_(v2wfpsf.kappa_),
|
||||
E_(v2wfpsf.E_),
|
||||
yPlusLam_(v2wfpsf.yPlusLam_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
fixedValueFvPatchField<scalar>(v2wfpsf, iF)
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
@ -183,6 +111,12 @@ void fWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
);
|
||||
const v2fBase& v2fModel = refCast<const v2fBase>(turbModel);
|
||||
|
||||
const nutWallFunctionFvPatchScalarField& nutw =
|
||||
refCast<const nutWallFunctionFvPatchScalarField>
|
||||
(
|
||||
turbModel.nut()().boundaryField()[patchi]
|
||||
);
|
||||
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
|
||||
const tmp<volScalarField> tk = turbModel.k();
|
||||
@ -197,7 +131,7 @@ void fWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
|
||||
const scalar Cmu25 = pow025(Cmu_);
|
||||
const scalar Cmu25 = pow025(nutw.Cmu());
|
||||
|
||||
scalarField& f = *this;
|
||||
|
||||
@ -210,7 +144,7 @@ void fWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
|
||||
scalar yPlus = uTau*y[facei]/nuw[facei];
|
||||
|
||||
if (yPlus > yPlusLam_)
|
||||
if (yPlus > nutw.yPlusLam())
|
||||
{
|
||||
scalar N = 6.0;
|
||||
scalar v2c = v2[celli];
|
||||
@ -232,22 +166,6 @@ void fWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
}
|
||||
|
||||
|
||||
void fWallFunctionFvPatchScalarField::evaluate
|
||||
(
|
||||
const Pstream::commsTypes commsType
|
||||
)
|
||||
{
|
||||
fixedValueFvPatchField<scalar>::evaluate(commsType);
|
||||
}
|
||||
|
||||
|
||||
void fWallFunctionFvPatchScalarField::write(Ostream& os) const
|
||||
{
|
||||
writeLocalEntries(os);
|
||||
fixedValueFvPatchField<scalar>::write(os);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
makePatchTypeField
|
||||
|
||||
@ -29,16 +29,10 @@ Description
|
||||
function condition for low- and high Reynolds number, turbulent flow cases
|
||||
|
||||
The model operates in two modes, based on the computed laminar-to-turbulent
|
||||
switch-over y+ value derived from kappa and E.
|
||||
switch-over y+ value derived from kappa and E specified in the corresponding
|
||||
nutWallFunction.
|
||||
|
||||
Usage
|
||||
\table
|
||||
Property | Description | Required | Default value
|
||||
Cmu | model coefficient | no | 0.09
|
||||
kappa | Von Karman constant | no | 0.41
|
||||
E | model coefficient | no | 9.8
|
||||
\endtable
|
||||
|
||||
Example of the boundary condition specification:
|
||||
\verbatim
|
||||
<patchName>
|
||||
@ -75,34 +69,6 @@ class fWallFunctionFvPatchScalarField
|
||||
:
|
||||
public fixedValueFvPatchField<scalar>
|
||||
{
|
||||
protected:
|
||||
|
||||
// Protected data
|
||||
|
||||
//- Cmu coefficient
|
||||
scalar Cmu_;
|
||||
|
||||
//- Von Karman constant
|
||||
scalar kappa_;
|
||||
|
||||
//- E coefficient
|
||||
scalar E_;
|
||||
|
||||
//- Y+ at the edge of the laminar sublayer
|
||||
scalar yPlusLam_;
|
||||
|
||||
|
||||
// Protected Member Functions
|
||||
|
||||
//- Check the type of the patch
|
||||
virtual void checkType();
|
||||
|
||||
//- Write local wall function variables
|
||||
virtual void writeLocalEntries(Ostream&) const;
|
||||
|
||||
//- Calculate the Y+ at the edge of the laminar sublayer
|
||||
scalar yPlusLam(const scalar kappa, const scalar E);
|
||||
|
||||
|
||||
public:
|
||||
|
||||
@ -178,15 +144,6 @@ public:
|
||||
|
||||
//- Update the coefficients associated with the patch field
|
||||
virtual void updateCoeffs();
|
||||
|
||||
//- Evaluate the patchField
|
||||
virtual void evaluate(const Pstream::commsTypes);
|
||||
|
||||
|
||||
// I-O
|
||||
|
||||
//- Write
|
||||
virtual void write(Ostream&) const;
|
||||
};
|
||||
|
||||
|
||||
|
||||
@ -24,50 +24,15 @@ License
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "kLowReWallFunctionFvPatchScalarField.H"
|
||||
#include "nutWallFunctionFvPatchScalarField.H"
|
||||
#include "turbulenceModel.H"
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "volFields.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
#include "wallFvPatch.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
|
||||
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
|
||||
|
||||
void kLowReWallFunctionFvPatchScalarField::checkType()
|
||||
{
|
||||
if (!isA<wallFvPatch>(patch()))
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Invalid wall function specification" << nl
|
||||
<< " Patch type for patch " << patch().name()
|
||||
<< " must be wall" << nl
|
||||
<< " Current patch type is " << patch().type() << nl << endl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
scalar kLowReWallFunctionFvPatchScalarField::yPlusLam
|
||||
(
|
||||
const scalar kappa,
|
||||
const scalar E
|
||||
)
|
||||
{
|
||||
scalar ypl = 11.0;
|
||||
|
||||
for (int i=0; i<10; i++)
|
||||
{
|
||||
ypl = log(max(E*ypl, 1))/kappa;
|
||||
}
|
||||
|
||||
return ypl;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
|
||||
@ -77,14 +42,8 @@ kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(p, iF),
|
||||
Cmu_(0.09),
|
||||
kappa_(0.41),
|
||||
E_(9.8),
|
||||
Ceps2_(1.9),
|
||||
yPlusLam_(yPlusLam(kappa_, E_))
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Ceps2_(1.9)
|
||||
{}
|
||||
|
||||
|
||||
kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
|
||||
@ -96,14 +55,8 @@ kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
|
||||
Cmu_(ptf.Cmu_),
|
||||
kappa_(ptf.kappa_),
|
||||
E_(ptf.E_),
|
||||
Ceps2_(ptf.Ceps2_),
|
||||
yPlusLam_(ptf.yPlusLam_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Ceps2_(ptf.Ceps2_)
|
||||
{}
|
||||
|
||||
|
||||
kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
|
||||
@ -114,14 +67,8 @@ kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(p, iF, dict),
|
||||
Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
|
||||
kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
|
||||
E_(dict.lookupOrDefault<scalar>("E", 9.8)),
|
||||
Ceps2_(dict.lookupOrDefault<scalar>("Ceps2", 1.9)),
|
||||
yPlusLam_(yPlusLam(kappa_, E_))
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Ceps2_(dict.lookupOrDefault<scalar>("Ceps2", 1.9))
|
||||
{}
|
||||
|
||||
|
||||
kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
|
||||
@ -130,14 +77,8 @@ kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(kwfpsf),
|
||||
Cmu_(kwfpsf.Cmu_),
|
||||
kappa_(kwfpsf.kappa_),
|
||||
E_(kwfpsf.E_),
|
||||
Ceps2_(kwfpsf.Ceps2_),
|
||||
yPlusLam_(kwfpsf.yPlusLam_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Ceps2_(kwfpsf.Ceps2_)
|
||||
{}
|
||||
|
||||
|
||||
kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
|
||||
@ -147,14 +88,8 @@ kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(kwfpsf, iF),
|
||||
Cmu_(kwfpsf.Cmu_),
|
||||
kappa_(kwfpsf.kappa_),
|
||||
E_(kwfpsf.E_),
|
||||
Ceps2_(kwfpsf.Ceps2_),
|
||||
yPlusLam_(kwfpsf.yPlusLam_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
Ceps2_(kwfpsf.Ceps2_)
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
@ -176,6 +111,13 @@ void kLowReWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
internalField().group()
|
||||
)
|
||||
);
|
||||
|
||||
const nutWallFunctionFvPatchScalarField& nutw =
|
||||
refCast<const nutWallFunctionFvPatchScalarField>
|
||||
(
|
||||
turbModel.nut()().boundaryField()[patchi]
|
||||
);
|
||||
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
|
||||
const tmp<volScalarField> tk = turbModel.k();
|
||||
@ -184,7 +126,7 @@ void kLowReWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
|
||||
const scalar Cmu25 = pow025(Cmu_);
|
||||
const scalar Cmu25 = pow025(nutw.Cmu());
|
||||
|
||||
scalarField& kw = *this;
|
||||
|
||||
@ -197,11 +139,11 @@ void kLowReWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
|
||||
scalar yPlus = uTau*y[facei]/nuw[facei];
|
||||
|
||||
if (yPlus > yPlusLam_)
|
||||
if (yPlus > nutw.yPlusLam())
|
||||
{
|
||||
scalar Ck = -0.416;
|
||||
scalar Bk = 8.366;
|
||||
kw[facei] = Ck/kappa_*log(yPlus) + Bk;
|
||||
kw[facei] = Ck/nutw.kappa()*log(yPlus) + Bk;
|
||||
}
|
||||
else
|
||||
{
|
||||
@ -222,20 +164,8 @@ void kLowReWallFunctionFvPatchScalarField::updateCoeffs()
|
||||
}
|
||||
|
||||
|
||||
void kLowReWallFunctionFvPatchScalarField::evaluate
|
||||
(
|
||||
const Pstream::commsTypes commsType
|
||||
)
|
||||
{
|
||||
fixedValueFvPatchField<scalar>::evaluate(commsType);
|
||||
}
|
||||
|
||||
|
||||
void kLowReWallFunctionFvPatchScalarField::write(Ostream& os) const
|
||||
{
|
||||
writeEntry(os, "Cmu", Cmu_);
|
||||
writeEntry(os, "kappa", kappa_);
|
||||
writeEntry(os, "E", E_);
|
||||
writeEntry(os, "Ceps2", Ceps2_);
|
||||
fixedValueFvPatchField<scalar>::write(os);
|
||||
}
|
||||
|
||||
@ -29,14 +29,12 @@ Description
|
||||
condition for low- and high-Reynolds number turbulent flow cases.
|
||||
|
||||
The model operates in two modes, based on the computed laminar-to-turbulent
|
||||
switch-over y+ value derived from kappa and E.
|
||||
switch-over y+ value derived from kappa and E specified in the corresponding
|
||||
nutWallFunction.
|
||||
|
||||
Usage
|
||||
\table
|
||||
Property | Description | Required | Default value
|
||||
Cmu | model coefficient | no | 0.09
|
||||
kappa | Von Karman constant | no | 0.41
|
||||
E | model coefficient | no | 9.8
|
||||
Ceps2 | model coefficient | no | 1.9
|
||||
\endtable
|
||||
|
||||
@ -78,30 +76,9 @@ protected:
|
||||
|
||||
// Protected data
|
||||
|
||||
//- Cmu coefficient
|
||||
scalar Cmu_;
|
||||
|
||||
//- Von Karman constant
|
||||
scalar kappa_;
|
||||
|
||||
//- E coefficient
|
||||
scalar E_;
|
||||
|
||||
//- Ceps2 coefficient
|
||||
scalar Ceps2_;
|
||||
|
||||
//- Y+ at the edge of the laminar sublayer
|
||||
scalar yPlusLam_;
|
||||
|
||||
|
||||
// Protected Member Functions
|
||||
|
||||
//- Check the type of the patch
|
||||
virtual void checkType();
|
||||
|
||||
//- Calculate the Y+ at the edge of the laminar sublayer
|
||||
scalar yPlusLam(const scalar kappa, const scalar E);
|
||||
|
||||
|
||||
public:
|
||||
|
||||
@ -178,9 +155,6 @@ public:
|
||||
//- Update the coefficients associated with the patch field
|
||||
virtual void updateCoeffs();
|
||||
|
||||
//- Evaluate the patchField
|
||||
virtual void evaluate(const Pstream::commsTypes);
|
||||
|
||||
|
||||
// I-O
|
||||
|
||||
|
||||
@ -24,27 +24,7 @@ License
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "kqRWallFunctionFvPatchField.H"
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
#include "wallFvPatch.H"
|
||||
|
||||
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||
|
||||
template<class Type>
|
||||
void Foam::kqRWallFunctionFvPatchField<Type>::checkType()
|
||||
{
|
||||
if (!isA<wallFvPatch>(this->patch()))
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Invalid wall function specification" << nl
|
||||
<< " Patch type for patch " << this->patch().name()
|
||||
<< " must be wall" << nl
|
||||
<< " Current patch type is " << this->patch().type()
|
||||
<< nl << endl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
@ -56,9 +36,7 @@ Foam::kqRWallFunctionFvPatchField<Type>::kqRWallFunctionFvPatchField
|
||||
)
|
||||
:
|
||||
zeroGradientFvPatchField<Type>(p, iF)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
template<class Type>
|
||||
@ -70,9 +48,7 @@ Foam::kqRWallFunctionFvPatchField<Type>::kqRWallFunctionFvPatchField
|
||||
)
|
||||
:
|
||||
zeroGradientFvPatchField<Type>(p, iF, dict)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
template<class Type>
|
||||
@ -85,9 +61,7 @@ Foam::kqRWallFunctionFvPatchField<Type>::kqRWallFunctionFvPatchField
|
||||
)
|
||||
:
|
||||
zeroGradientFvPatchField<Type>(ptf, p, iF, mapper)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
template<class Type>
|
||||
@ -97,9 +71,7 @@ Foam::kqRWallFunctionFvPatchField<Type>::kqRWallFunctionFvPatchField
|
||||
)
|
||||
:
|
||||
zeroGradientFvPatchField<Type>(tkqrwfpf)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
template<class Type>
|
||||
@ -110,9 +82,7 @@ Foam::kqRWallFunctionFvPatchField<Type>::kqRWallFunctionFvPatchField
|
||||
)
|
||||
:
|
||||
zeroGradientFvPatchField<Type>(tkqrwfpf, iF)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
@ -68,12 +68,6 @@ class kqRWallFunctionFvPatchField
|
||||
public zeroGradientFvPatchField<Type>
|
||||
{
|
||||
|
||||
// Private Member Functions
|
||||
|
||||
//- Check the type of the patch
|
||||
void checkType();
|
||||
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
|
||||
@ -35,6 +35,7 @@ Usage
|
||||
<patchName>
|
||||
{
|
||||
type nutLowReWallFunction;
|
||||
value uniform 0;
|
||||
}
|
||||
\endverbatim
|
||||
|
||||
|
||||
@ -44,6 +44,7 @@ Usage
|
||||
roughnessHeight 1e-5;
|
||||
roughnessConstant 0.5;
|
||||
roughnessFactor 1;
|
||||
value uniform 0;
|
||||
}
|
||||
\endverbatim
|
||||
|
||||
|
||||
@ -48,6 +48,7 @@ Usage
|
||||
<patchName>
|
||||
{
|
||||
type nutUSpaldingWallFunction;
|
||||
value uniform 0;
|
||||
}
|
||||
\endverbatim
|
||||
|
||||
|
||||
@ -42,6 +42,7 @@ Usage
|
||||
{
|
||||
type nutTabulatedWallFunction;
|
||||
uPlusTable myUPlusTable;
|
||||
value uniform 0;
|
||||
}
|
||||
\endverbatim
|
||||
|
||||
|
||||
@ -34,6 +34,7 @@ Usage
|
||||
<patchName>
|
||||
{
|
||||
type nutUWallFunction;
|
||||
value uniform 0;
|
||||
}
|
||||
\endverbatim
|
||||
|
||||
|
||||
@ -42,7 +42,7 @@ Usage
|
||||
<patchName>
|
||||
{
|
||||
type nutWallFunction;
|
||||
value uniform 0.0;
|
||||
value uniform 0;
|
||||
}
|
||||
\endverbatim
|
||||
|
||||
@ -158,6 +158,24 @@ public:
|
||||
|
||||
// Member functions
|
||||
|
||||
//- Return Cmu
|
||||
scalar Cmu() const
|
||||
{
|
||||
return Cmu_;
|
||||
}
|
||||
|
||||
//- Return kappa
|
||||
scalar kappa() const
|
||||
{
|
||||
return kappa_;
|
||||
}
|
||||
|
||||
//- Return E
|
||||
scalar E() const
|
||||
{
|
||||
return E_;
|
||||
}
|
||||
|
||||
//- Calculate the Y+ at the edge of the laminar sublayer
|
||||
static scalar yPlusLam(const scalar kappa, const scalar E);
|
||||
|
||||
|
||||
@ -48,6 +48,7 @@ Usage
|
||||
type nutkRoughWallFunction;
|
||||
Ks uniform 0;
|
||||
Cs uniform 0.5;
|
||||
value uniform 0;
|
||||
}
|
||||
\endverbatim
|
||||
|
||||
|
||||
@ -34,6 +34,7 @@ Usage
|
||||
<patchName>
|
||||
{
|
||||
type nutkWallFunction;
|
||||
value uniform 0;
|
||||
}
|
||||
\endverbatim
|
||||
|
||||
|
||||
@ -26,10 +26,7 @@ License
|
||||
#include "omegaWallFunctionFvPatchScalarField.H"
|
||||
#include "nutWallFunctionFvPatchScalarField.H"
|
||||
#include "turbulenceModel.H"
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "fvMatrix.H"
|
||||
#include "volFields.H"
|
||||
#include "wallFvPatch.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
@ -43,30 +40,6 @@ scalar omegaWallFunctionFvPatchScalarField::tolerance_ = 1e-5;
|
||||
|
||||
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
|
||||
|
||||
void omegaWallFunctionFvPatchScalarField::checkType()
|
||||
{
|
||||
if (!isA<wallFvPatch>(patch()))
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Invalid wall function specification" << nl
|
||||
<< " Patch type for patch " << patch().name()
|
||||
<< " must be wall" << nl
|
||||
<< " Current patch type is " << patch().type() << nl << endl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
void omegaWallFunctionFvPatchScalarField::writeLocalEntries(Ostream& os) const
|
||||
{
|
||||
writeEntry(os, "Cmu", Cmu_);
|
||||
writeEntry(os, "kappa", kappa_);
|
||||
writeEntry(os, "E", E_);
|
||||
writeEntry(os, "beta1", beta1_);
|
||||
writeEntry(os, "blended", blended_);
|
||||
}
|
||||
|
||||
|
||||
void omegaWallFunctionFvPatchScalarField::setMaster()
|
||||
{
|
||||
if (master_ != -1)
|
||||
@ -216,10 +189,13 @@ void omegaWallFunctionFvPatchScalarField::calculate
|
||||
{
|
||||
const label patchi = patch.index();
|
||||
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
const nutWallFunctionFvPatchScalarField& nutw =
|
||||
refCast<const nutWallFunctionFvPatchScalarField>
|
||||
(
|
||||
turbModel.nut()().boundaryField()[patchi]
|
||||
);
|
||||
|
||||
const scalar Cmu25 = pow025(Cmu_);
|
||||
const scalar Cmu5 = sqrt(Cmu_);
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
|
||||
const tmp<volScalarField> tk = turbModel.k();
|
||||
const volScalarField& k = tk();
|
||||
@ -227,9 +203,6 @@ void omegaWallFunctionFvPatchScalarField::calculate
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
|
||||
const tmp<scalarField> tnutw = turbModel.nut(patchi);
|
||||
const scalarField& nutw = tnutw();
|
||||
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
|
||||
const scalarField magGradUw(mag(Uw.snGrad()));
|
||||
@ -237,6 +210,9 @@ void omegaWallFunctionFvPatchScalarField::calculate
|
||||
const FieldType& G =
|
||||
db().lookupObject<FieldType>(turbModel.GName());
|
||||
|
||||
const scalar Cmu25 = pow025(nutw.Cmu());
|
||||
const scalar Cmu5 = sqrt(nutw.Cmu());
|
||||
|
||||
// Set omega and G
|
||||
forAll(nutw, facei)
|
||||
{
|
||||
@ -245,7 +221,7 @@ void omegaWallFunctionFvPatchScalarField::calculate
|
||||
|
||||
const scalar Rey = y[facei]*sqrt(k[celli])/nuw[facei];
|
||||
const scalar yPlus = Cmu25*Rey;
|
||||
const scalar uPlus = (1/kappa_)*log(E_*yPlus);
|
||||
const scalar uPlus = (1/nutw.kappa())*log(nutw.E()*yPlus);
|
||||
|
||||
if (blended_)
|
||||
{
|
||||
@ -258,7 +234,7 @@ void omegaWallFunctionFvPatchScalarField::calculate
|
||||
);
|
||||
|
||||
const scalar omegaVis = 6*nuw[facei]/(beta1_*sqr(y[facei]));
|
||||
const scalar omegaLog = uStar/(Cmu5*kappa_*y[facei]);
|
||||
const scalar omegaLog = uStar/(Cmu5*nutw.kappa()*y[facei]);
|
||||
|
||||
omega0[celli] += w*(lamFrac*omegaVis + turbFrac*omegaLog);
|
||||
|
||||
@ -269,12 +245,12 @@ void omegaWallFunctionFvPatchScalarField::calculate
|
||||
|
||||
+ turbFrac
|
||||
*sqr(uStar*magGradUw[facei]*y[facei]/uPlus)
|
||||
/(nuw[facei]*kappa_*yPlus)
|
||||
/(nuw[facei]*nutw.kappa()*yPlus)
|
||||
);
|
||||
}
|
||||
else
|
||||
{
|
||||
if (yPlus < yPlusLam_)
|
||||
if (yPlus < nutw.yPlusLam())
|
||||
{
|
||||
const scalar omegaVis = 6*nuw[facei]/(beta1_*sqr(y[facei]));
|
||||
|
||||
@ -285,14 +261,14 @@ void omegaWallFunctionFvPatchScalarField::calculate
|
||||
else
|
||||
{
|
||||
const scalar uStar = sqrt(Cmu5*k[celli]);
|
||||
const scalar omegaLog = uStar/(Cmu5*kappa_*y[facei]);
|
||||
const scalar omegaLog = uStar/(Cmu5*nutw.kappa()*y[facei]);
|
||||
|
||||
omega0[celli] += w*omegaLog;
|
||||
|
||||
G0[celli] +=
|
||||
w*
|
||||
sqr(uStar*magGradUw[facei]*y[facei]/uPlus)
|
||||
/(nuw[facei]*kappa_*yPlus);
|
||||
/(nuw[facei]*nutw.kappa()*yPlus);
|
||||
}
|
||||
}
|
||||
}
|
||||
@ -308,20 +284,14 @@ omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(p, iF),
|
||||
Cmu_(0.09),
|
||||
kappa_(0.41),
|
||||
E_(9.8),
|
||||
beta1_(0.075),
|
||||
blended_(false),
|
||||
yPlusLam_(nutWallFunctionFvPatchScalarField::yPlusLam(kappa_, E_)),
|
||||
G_(),
|
||||
omega_(),
|
||||
initialised_(false),
|
||||
master_(-1),
|
||||
cornerWeights_()
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
@ -332,20 +302,14 @@ omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(p, iF, dict),
|
||||
Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
|
||||
kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
|
||||
E_(dict.lookupOrDefault<scalar>("E", 9.8)),
|
||||
beta1_(dict.lookupOrDefault<scalar>("beta1", 0.075)),
|
||||
blended_(dict.lookupOrDefault<Switch>("blended", false)),
|
||||
yPlusLam_(nutWallFunctionFvPatchScalarField::yPlusLam(kappa_, E_)),
|
||||
G_(),
|
||||
omega_(),
|
||||
initialised_(false),
|
||||
master_(-1),
|
||||
cornerWeights_()
|
||||
{
|
||||
checkType();
|
||||
|
||||
// apply zero-gradient condition on start-up
|
||||
this->operator==(patchInternalField());
|
||||
}
|
||||
@ -360,20 +324,14 @@ omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
|
||||
Cmu_(ptf.Cmu_),
|
||||
kappa_(ptf.kappa_),
|
||||
E_(ptf.E_),
|
||||
beta1_(ptf.beta1_),
|
||||
blended_(ptf.blended_),
|
||||
yPlusLam_(ptf.yPlusLam_),
|
||||
G_(),
|
||||
omega_(),
|
||||
initialised_(false),
|
||||
master_(-1),
|
||||
cornerWeights_()
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
@ -382,20 +340,14 @@ omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(owfpsf),
|
||||
Cmu_(owfpsf.Cmu_),
|
||||
kappa_(owfpsf.kappa_),
|
||||
E_(owfpsf.E_),
|
||||
beta1_(owfpsf.beta1_),
|
||||
blended_(owfpsf.blended_),
|
||||
yPlusLam_(owfpsf.yPlusLam_),
|
||||
G_(),
|
||||
omega_(),
|
||||
initialised_(false),
|
||||
master_(-1),
|
||||
cornerWeights_()
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
@ -405,20 +357,14 @@ omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(owfpsf, iF),
|
||||
Cmu_(owfpsf.Cmu_),
|
||||
kappa_(owfpsf.kappa_),
|
||||
E_(owfpsf.E_),
|
||||
beta1_(owfpsf.beta1_),
|
||||
blended_(owfpsf.blended_),
|
||||
yPlusLam_(owfpsf.yPlusLam_),
|
||||
G_(),
|
||||
omega_(),
|
||||
initialised_(false),
|
||||
master_(-1),
|
||||
cornerWeights_()
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
@ -632,7 +578,8 @@ void omegaWallFunctionFvPatchScalarField::manipulateMatrix
|
||||
|
||||
void omegaWallFunctionFvPatchScalarField::write(Ostream& os) const
|
||||
{
|
||||
writeLocalEntries(os);
|
||||
writeEntry(os, "beta1", beta1_);
|
||||
writeEntry(os, "blended", blended_);
|
||||
fixedValueFvPatchField<scalar>::write(os);
|
||||
}
|
||||
|
||||
|
||||
@ -51,18 +51,16 @@ Description
|
||||
\endverbatim
|
||||
|
||||
or switched between these values based on the laminar-to-turbulent y+ value
|
||||
derived from kappa and E. Recent tests have shown that the standard
|
||||
switching method provides more accurate results for 10 < y+ < 30 when used
|
||||
with high Reynolds number wall-functions and both methods provide accurate
|
||||
results when used with continuous wall-functions. Based on this the
|
||||
standard switching method is used by default.
|
||||
derived from kappa and E specified in the corresponding nutWallFunction.
|
||||
Recent tests have shown that the standard switching method provides more
|
||||
accurate results for 10 < y+ < 30 when used with high Reynolds number
|
||||
wall-functions and both methods provide accurate results when used with
|
||||
continuous wall-functions. Based on this the standard switching method is
|
||||
used by default.
|
||||
|
||||
Usage
|
||||
\table
|
||||
Property | Description | Required | Default value
|
||||
Cmu | Model coefficient | no | 0.09
|
||||
kappa | von Karman constant | no | 0.41
|
||||
E | Model coefficient | no | 9.8
|
||||
beta1 | Model coefficient | no | 0.075
|
||||
blended | Blending switch | no | false
|
||||
\endtable
|
||||
@ -111,24 +109,12 @@ protected:
|
||||
//- Tolerance used in weighted calculations
|
||||
static scalar tolerance_;
|
||||
|
||||
//- Cmu coefficient
|
||||
scalar Cmu_;
|
||||
|
||||
//- Von Karman constant
|
||||
scalar kappa_;
|
||||
|
||||
//- E coefficient
|
||||
scalar E_;
|
||||
|
||||
//- beta1 coefficient
|
||||
scalar beta1_;
|
||||
|
||||
//- Blending switch (defaults to false)
|
||||
Switch blended_;
|
||||
|
||||
//- y+ at the edge of the laminar sublayer
|
||||
scalar yPlusLam_;
|
||||
|
||||
//- Local copy of turbulence G field
|
||||
scalarField G_;
|
||||
|
||||
@ -147,12 +133,6 @@ protected:
|
||||
|
||||
// Protected Member Functions
|
||||
|
||||
//- Check the type of the patch
|
||||
virtual void checkType();
|
||||
|
||||
//- Write local wall function variables
|
||||
virtual void writeLocalEntries(Ostream&) const;
|
||||
|
||||
//- Set the master patch - master is responsible for updating all
|
||||
// wall function patches
|
||||
virtual void setMaster();
|
||||
|
||||
@ -24,11 +24,9 @@ License
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "v2WallFunctionFvPatchScalarField.H"
|
||||
#include "nutWallFunctionFvPatchScalarField.H"
|
||||
#include "turbulenceModel.H"
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "volFields.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
#include "wallFvPatch.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -37,47 +35,6 @@ namespace Foam
|
||||
namespace RASModels
|
||||
{
|
||||
|
||||
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
|
||||
|
||||
void v2WallFunctionFvPatchScalarField::checkType()
|
||||
{
|
||||
if (!isA<wallFvPatch>(patch()))
|
||||
{
|
||||
FatalErrorInFunction
|
||||
<< "Invalid wall function specification" << nl
|
||||
<< " Patch type for patch " << patch().name()
|
||||
<< " must be wall" << nl
|
||||
<< " Current patch type is " << patch().type() << nl << endl
|
||||
<< abort(FatalError);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
void v2WallFunctionFvPatchScalarField::writeLocalEntries(Ostream& os) const
|
||||
{
|
||||
writeEntry(os, "Cmu", Cmu_);
|
||||
writeEntry(os, "kappa", kappa_);
|
||||
writeEntry(os, "E", E_);
|
||||
}
|
||||
|
||||
|
||||
scalar v2WallFunctionFvPatchScalarField::yPlusLam
|
||||
(
|
||||
const scalar kappa,
|
||||
const scalar E
|
||||
)
|
||||
{
|
||||
scalar ypl = 11.0;
|
||||
|
||||
for (int i=0; i<10; i++)
|
||||
{
|
||||
ypl = log(max(E*ypl, 1))/kappa;
|
||||
}
|
||||
|
||||
return ypl;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
v2WallFunctionFvPatchScalarField::v2WallFunctionFvPatchScalarField
|
||||
@ -86,14 +43,8 @@ v2WallFunctionFvPatchScalarField::v2WallFunctionFvPatchScalarField
|
||||
const DimensionedField<scalar, volMesh>& iF
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(p, iF),
|
||||
Cmu_(0.09),
|
||||
kappa_(0.41),
|
||||
E_(9.8),
|
||||
yPlusLam_(yPlusLam(kappa_, E_))
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
fixedValueFvPatchField<scalar>(p, iF)
|
||||
{}
|
||||
|
||||
|
||||
v2WallFunctionFvPatchScalarField::v2WallFunctionFvPatchScalarField
|
||||
@ -104,14 +55,8 @@ v2WallFunctionFvPatchScalarField::v2WallFunctionFvPatchScalarField
|
||||
const fvPatchFieldMapper& mapper
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
|
||||
Cmu_(ptf.Cmu_),
|
||||
kappa_(ptf.kappa_),
|
||||
E_(ptf.E_),
|
||||
yPlusLam_(ptf.yPlusLam_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
fixedValueFvPatchField<scalar>(ptf, p, iF, mapper)
|
||||
{}
|
||||
|
||||
|
||||
v2WallFunctionFvPatchScalarField::v2WallFunctionFvPatchScalarField
|
||||
@ -121,14 +66,8 @@ v2WallFunctionFvPatchScalarField::v2WallFunctionFvPatchScalarField
|
||||
const dictionary& dict
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(p, iF, dict),
|
||||
Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
|
||||
kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
|
||||
E_(dict.lookupOrDefault<scalar>("E", 9.8)),
|
||||
yPlusLam_(yPlusLam(kappa_, E_))
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
fixedValueFvPatchField<scalar>(p, iF, dict)
|
||||
{}
|
||||
|
||||
|
||||
v2WallFunctionFvPatchScalarField::v2WallFunctionFvPatchScalarField
|
||||
@ -136,14 +75,8 @@ v2WallFunctionFvPatchScalarField::v2WallFunctionFvPatchScalarField
|
||||
const v2WallFunctionFvPatchScalarField& v2wfpsf
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(v2wfpsf),
|
||||
Cmu_(v2wfpsf.Cmu_),
|
||||
kappa_(v2wfpsf.kappa_),
|
||||
E_(v2wfpsf.E_),
|
||||
yPlusLam_(v2wfpsf.yPlusLam_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
fixedValueFvPatchField<scalar>(v2wfpsf)
|
||||
{}
|
||||
|
||||
|
||||
v2WallFunctionFvPatchScalarField::v2WallFunctionFvPatchScalarField
|
||||
@ -152,14 +85,8 @@ v2WallFunctionFvPatchScalarField::v2WallFunctionFvPatchScalarField
|
||||
const DimensionedField<scalar, volMesh>& iF
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchField<scalar>(v2wfpsf, iF),
|
||||
Cmu_(v2wfpsf.Cmu_),
|
||||
kappa_(v2wfpsf.kappa_),
|
||||
E_(v2wfpsf.E_),
|
||||
yPlusLam_(v2wfpsf.yPlusLam_)
|
||||
{
|
||||
checkType();
|
||||
}
|
||||
fixedValueFvPatchField<scalar>(v2wfpsf, iF)
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
@ -181,6 +108,13 @@ void v2WallFunctionFvPatchScalarField::updateCoeffs()
|
||||
internalField().group()
|
||||
)
|
||||
);
|
||||
|
||||
const nutWallFunctionFvPatchScalarField& nutw =
|
||||
refCast<const nutWallFunctionFvPatchScalarField>
|
||||
(
|
||||
turbModel.nut()().boundaryField()[patchi]
|
||||
);
|
||||
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
|
||||
const tmp<volScalarField> tk = turbModel.k();
|
||||
@ -189,7 +123,7 @@ void v2WallFunctionFvPatchScalarField::updateCoeffs()
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
|
||||
const scalar Cmu25 = pow025(Cmu_);
|
||||
const scalar Cmu25 = pow025(nutw.Cmu());
|
||||
|
||||
scalarField& v2 = *this;
|
||||
|
||||
@ -202,11 +136,11 @@ void v2WallFunctionFvPatchScalarField::updateCoeffs()
|
||||
|
||||
scalar yPlus = uTau*y[facei]/nuw[facei];
|
||||
|
||||
if (yPlus > yPlusLam_)
|
||||
if (yPlus > nutw.yPlusLam())
|
||||
{
|
||||
scalar Cv2 = 0.193;
|
||||
scalar Bv2 = -0.94;
|
||||
v2[facei] = Cv2/kappa_*log(yPlus) + Bv2;
|
||||
v2[facei] = Cv2/nutw.kappa()*log(yPlus) + Bv2;
|
||||
}
|
||||
else
|
||||
{
|
||||
@ -223,22 +157,6 @@ void v2WallFunctionFvPatchScalarField::updateCoeffs()
|
||||
}
|
||||
|
||||
|
||||
void v2WallFunctionFvPatchScalarField::evaluate
|
||||
(
|
||||
const Pstream::commsTypes commsType
|
||||
)
|
||||
{
|
||||
fixedValueFvPatchField<scalar>::evaluate(commsType);
|
||||
}
|
||||
|
||||
|
||||
void v2WallFunctionFvPatchScalarField::write(Ostream& os) const
|
||||
{
|
||||
writeLocalEntries(os);
|
||||
fixedValueFvPatchField<scalar>::write(os);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
makePatchTypeField
|
||||
|
||||
@ -30,17 +30,11 @@ Description
|
||||
cases.
|
||||
|
||||
The model operates in two modes, based on the computed laminar-to-turbulent
|
||||
switch-over y+ value derived from kappa and E.
|
||||
switch-over y+ value derived from kappa and E specified in the corresponding
|
||||
nutWallFunction.
|
||||
|
||||
|
||||
Usage
|
||||
\table
|
||||
Property | Description | Required | Default value
|
||||
Cmu | model coefficient | no | 0.09
|
||||
kappa | Von Karman constant | no | 0.41
|
||||
E | model coefficient | no | 9.8
|
||||
\endtable
|
||||
|
||||
Example of the boundary condition specification:
|
||||
\verbatim
|
||||
<patchName>
|
||||
@ -77,34 +71,6 @@ class v2WallFunctionFvPatchScalarField
|
||||
:
|
||||
public fixedValueFvPatchField<scalar>
|
||||
{
|
||||
protected:
|
||||
|
||||
// Protected data
|
||||
|
||||
//- Cmu coefficient
|
||||
scalar Cmu_;
|
||||
|
||||
//- Von Karman constant
|
||||
scalar kappa_;
|
||||
|
||||
//- E coefficient
|
||||
scalar E_;
|
||||
|
||||
//- Y+ at the edge of the laminar sublayer
|
||||
scalar yPlusLam_;
|
||||
|
||||
|
||||
// Protected Member Functions
|
||||
|
||||
//- Check the type of the patch
|
||||
virtual void checkType();
|
||||
|
||||
//- Write local wall function variables
|
||||
virtual void writeLocalEntries(Ostream&) const;
|
||||
|
||||
//- Calculate the Y+ at the edge of the laminar sublayer
|
||||
scalar yPlusLam(const scalar kappa, const scalar E);
|
||||
|
||||
|
||||
public:
|
||||
|
||||
@ -180,15 +146,6 @@ public:
|
||||
|
||||
//- Update the coefficients associated with the patch field
|
||||
virtual void updateCoeffs();
|
||||
|
||||
//- Evaluate the patchField
|
||||
virtual void evaluate(const Pstream::commsTypes);
|
||||
|
||||
|
||||
// I-O
|
||||
|
||||
//- Write
|
||||
virtual void write(Ostream&) const;
|
||||
};
|
||||
|
||||
|
||||
|
||||
@ -35,7 +35,7 @@ boundaryField
|
||||
|
||||
plate
|
||||
{
|
||||
type fixedValue;
|
||||
type nutkWallFunction;
|
||||
value uniform 0;
|
||||
}
|
||||
|
||||
|
||||
Reference in New Issue
Block a user