The k-epsilon-phit-f turbulence closure model for incompressible and compressible flows. More...
Public Types | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
![]() | |
typedef RASModel< BasicTurbulenceModel > ::alphaField | alphaField |
typedef RASModel< BasicTurbulenceModel > ::rhoField | rhoField |
typedef RASModel< BasicTurbulenceModel > ::transportModel | transportModel |
![]() | |
typedef RASModel< BasicTurbulenceModel > ::alphaField | alphaField |
typedef RASModel< BasicTurbulenceModel > ::rhoField | rhoField |
typedef RASModel< BasicTurbulenceModel > ::transportModel | transportModel |
![]() | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
TypeName ("kEpsilonPhitF") | |
Runtime type information. More... | |
kEpsilonPhitF (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName) | |
Construct from components. More... | |
virtual | ~kEpsilonPhitF ()=default |
Destructor. More... | |
virtual bool | read () |
Re-read model coefficients if they have changed. More... | |
tmp< volScalarField > | DkEff () const |
Return the effective diffusivity for k (LUU:Eq. 3) More... | |
tmp< volScalarField > | DepsilonEff () const |
Return the effective diffusivity for epsilon (LUU:Eq. 4) More... | |
tmp< volScalarField > | DphitEff () const |
Return the effective diffusivity for phit (LUU:Eq. 17) More... | |
virtual tmp< volScalarField > | k () const |
Return the turbulent kinetic energy field. More... | |
virtual tmp< volScalarField > | epsilon () const |
Return the turbulent kinetic energy dissipation rate field. More... | |
virtual tmp< volScalarField > | phit () const |
Return the normalised wall-normal fluctuating velocity scale field. More... | |
virtual tmp< volScalarField > | f () const |
Return the elliptic relaxation factor field. More... | |
virtual void | correct () |
Solve the transport equations and correct the turbulent viscosity. More... | |
![]() | |
eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName) | |
Construct from components. More... | |
virtual | ~eddyViscosity ()=default |
Destructor. More... | |
virtual tmp< volScalarField > | nut () const |
Return the turbulence viscosity. More... | |
virtual tmp< scalarField > | nut (const label patchi) const |
Return the turbulence viscosity on patch. More... | |
virtual tmp< volSymmTensorField > | R () const |
Return the Reynolds stress tensor. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
![]() | |
linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName) | |
Construct from components. More... | |
virtual | ~linearViscousStress ()=default |
Destructor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff () const |
Return the effective stress tensor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
Return the effective stress tensor based on a given velocity field. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
Return the source term for the momentum equation. More... | |
![]() | |
TypeName ("RAS") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (autoPtr, RASModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)) | |
RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName) | |
Construct from components. More... | |
virtual | ~RASModel ()=default |
Destructor. More... | |
const dimensionedScalar & | kMin () const |
Return the lower allowable limit for k (default: SMALL) More... | |
const dimensionedScalar & | epsilonMin () const |
Return the lower allowable limit for epsilon (default: SMALL) More... | |
const dimensionedScalar & | omegaMin () const |
Return the lower allowable limit for omega (default: SMALL) More... | |
dimensionedScalar & | kMin () |
Allow kMin to be changed. More... | |
dimensionedScalar & | epsilonMin () |
Allow epsilonMin to be changed. More... | |
dimensionedScalar & | omegaMin () |
Allow omegaMin to be changed. More... | |
virtual const dictionary & | coeffDict () const |
Const access to the coefficients dictionary. More... | |
virtual tmp< volScalarField > | nuEff () const |
Return the effective viscosity. More... | |
virtual tmp< scalarField > | nuEff (const label patchi) const |
Return the effective viscosity on patch. More... | |
virtual tmp< volScalarField > | omega () const |
Return the specific dissipation rate. More... | |
Protected Member Functions | |
virtual void | correctNut () |
Update nut with the latest available k, phit, and T. More... | |
tmp< volScalarField > | Ts () const |
Return the turbulent time scale, T. More... | |
tmp< volScalarField > | Ls () const |
Return the turbulent length scale, L. More... | |
![]() | |
virtual void | printCoeffs (const word &type) |
Print model coefficients. More... | |
RASModel (const RASModel &)=delete | |
No copy construct. More... | |
void | operator= (const RASModel &)=delete |
No copy assignment. More... | |
Additional Inherited Members | |
![]() | |
static autoPtr< RASModel > | New (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName) |
Return a reference to the selected RAS model. More... | |
The k-epsilon-phit-f turbulence closure model for incompressible and compressible flows.
The model is a three-transport-equation linear-eddy-viscosity turbulence closure model alongside an elliptic relaxation equation.
Input fields
k | : | Turbulent kinetic energy [m2/s2] |
epsilon | : | Turbulent kinetic energy dissipation rate [m2/s3] |
phit | : | Normalised wall-normal fluctuating velocity scale [-] |
f | : | Elliptic relaxation factor [1/s] |
Reference:
Standard model (Tag:LUU): Laurence, D. R., Uribe, J. C., & Utyuzhnikov, S. V. (2005). A robust formulation of the v2-f model. Flow, Turbulence and Combustion, 73(3-4), 169–185. DOI:10.1007/s10494-005-1974-8
The default model coefficients are (LUU:Eqs. 19-20):
kEpsilonPhitFCoeffs { includeNu true; // include nu in (LUU: Eq. 17), see Notes Cmu 0.22; // Turbulent viscosity constant Ceps1a 1.4; // Model constant for epsilon Ceps1b 1.0; // Model constant for epsilon Ceps1c 0.05; // Model constant for epsilon Ceps2 1.9; // Model constant for epsilon Cf1 1.4; // Model constant for f Cf2 0.3; // Model constant for f CL 0.25; // Model constant for L Ceta 110.0; // Model constant for L CT 6.0; // Model constant for T sigmaK 1.0; // Turbulent Prandtl number for k sigmaEps 1.3; // Turbulent Prandtl number for epsilon sigmaPhit 1.0; // Turbulent Prandtl number for phit = sigmaK }
Including nu
in DphitEff
even though it is not present in (LUU:Eq. 17) provided higher level of resemblance to benchmarks for the tests considered, particularly for the peak skin friction (yet, pressure-related predictions were unaffected). Users can switch off nu
in DphitEff
by using includeNu
entry in kEpsilonPhitFCoeffs
as shown above in order to follow the reference paper thereat. includeNu
is left true
by default. See GitLab issue #1560.
Definition at line 131 of file kEpsilonPhitF.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 207 of file kEpsilonPhitF.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 208 of file kEpsilonPhitF.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 209 of file kEpsilonPhitF.H.
kEpsilonPhitF | ( | const alphaField & | alpha, |
const rhoField & | rho, | ||
const volVectorField & | U, | ||
const surfaceScalarField & | alphaRhoPhi, | ||
const surfaceScalarField & | phi, | ||
const transportModel & | transport, | ||
const word & | propertiesName = turbulenceModel::propertiesName , |
||
const word & | type = typeName |
||
) |
Construct from components.
Definition at line 100 of file kEpsilonPhitF.C.
References Foam::bound(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::mag(), Foam::nl, and Foam::type().
|
virtualdefault |
Destructor.
|
protectedvirtual |
Update nut with the latest available k, phit, and T.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 42 of file kEpsilonPhitF.C.
References optionList::correct(), and options::New().
|
protected |
Return the turbulent time scale, T.
Definition at line 54 of file kEpsilonPhitF.C.
References Foam::max(), nu, Foam::sqrt(), and Foam::Zero.
|
protected |
Return the turbulent length scale, L.
Definition at line 74 of file kEpsilonPhitF.C.
References Foam::max(), nu, Foam::pow(), Foam::pow025(), Foam::pow3(), and Foam::Zero.
TypeName | ( | "kEpsilonPhitF< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 348 of file kEpsilonPhitF.C.
References Foam::read().
|
inline |
Return the effective diffusivity for k (LUU:Eq. 3)
Definition at line 242 of file kEpsilonPhitF.H.
References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.
|
inline |
Return the effective diffusivity for epsilon (LUU:Eq. 4)
Definition at line 255 of file kEpsilonPhitF.H.
References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.
|
inline |
Return the effective diffusivity for phit (LUU:Eq. 17)
Definition at line 268 of file kEpsilonPhitF.H.
References kEpsilonPhitF< BasicTurbulenceModel >::includeNu_, tmp< T >::New(), nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.
|
inlinevirtual |
Return the turbulent kinetic energy field.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 282 of file kEpsilonPhitF.H.
References kEpsilonPhitF< BasicTurbulenceModel >::k_.
|
inlinevirtual |
Return the turbulent kinetic energy dissipation rate field.
Reimplemented from RASModel< BasicTurbulenceModel >.
Definition at line 288 of file kEpsilonPhitF.H.
References kEpsilonPhitF< BasicTurbulenceModel >::epsilon_.
|
inlinevirtual |
Return the normalised wall-normal fluctuating velocity scale field.
Definition at line 294 of file kEpsilonPhitF.H.
References kEpsilonPhitF< BasicTurbulenceModel >::phit_.
|
inlinevirtual |
Return the elliptic relaxation factor field.
Definition at line 300 of file kEpsilonPhitF.H.
References kEpsilonPhitF< BasicTurbulenceModel >::f_.
|
virtual |
Solve the transport equations and correct the turbulent viscosity.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 375 of file kEpsilonPhitF.C.
References Foam::fvc::absolute(), Foam::constant::atomic::alpha, Foam::bound(), tmp< T >::clear(), correct(), Foam::fvm::ddt(), Foam::dev(), Foam::fvm::div(), Foam::fvc::div(), divU, fvOptions, Foam::constant::universal::G, Foam::fvc::grad(), Foam::fvc::laplacian(), Foam::fvm::laplacian(), options::New(), nu, nut, phi, tmp< T >::ref(), rho, Foam::solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::sqrt(), Foam::fvm::SuSp(), Foam::symm(), Foam::type(), and U.
|
protected |
Definition at line 148 of file kEpsilonPhitF.H.
Referenced by kEpsilonPhitF< BasicTurbulenceModel >::DphitEff().
|
protected |
Definition at line 152 of file kEpsilonPhitF.H.
|
protected |
Definition at line 153 of file kEpsilonPhitF.H.
|
protected |
Definition at line 154 of file kEpsilonPhitF.H.
|
protected |
Definition at line 155 of file kEpsilonPhitF.H.
|
protected |
Definition at line 156 of file kEpsilonPhitF.H.
|
protected |
Definition at line 157 of file kEpsilonPhitF.H.
|
protected |
Definition at line 158 of file kEpsilonPhitF.H.
|
protected |
Definition at line 159 of file kEpsilonPhitF.H.
|
protected |
Definition at line 160 of file kEpsilonPhitF.H.
|
protected |
Definition at line 161 of file kEpsilonPhitF.H.
|
protected |
Definition at line 162 of file kEpsilonPhitF.H.
|
protected |
Definition at line 163 of file kEpsilonPhitF.H.
|
protected |
Definition at line 164 of file kEpsilonPhitF.H.
|
protected |
Turbulent kinetic energy [m2/s2].
Definition at line 170 of file kEpsilonPhitF.H.
Referenced by kEpsilonPhitF< BasicTurbulenceModel >::k().
|
protected |
Turbulent kinetic energy dissipation rate [m2/s3].
Definition at line 173 of file kEpsilonPhitF.H.
Referenced by kEpsilonPhitF< BasicTurbulenceModel >::epsilon().
|
protected |
Normalised wall-normal fluctuating velocity scale [-].
Definition at line 176 of file kEpsilonPhitF.H.
Referenced by kEpsilonPhitF< BasicTurbulenceModel >::phit().
|
protected |
Elliptic relaxation factor [1/s].
Definition at line 179 of file kEpsilonPhitF.H.
Referenced by kEpsilonPhitF< BasicTurbulenceModel >::f().
|
protected |
Turbulent time scale [s].
Definition at line 182 of file kEpsilonPhitF.H.
|
protected |
Definition at line 187 of file kEpsilonPhitF.H.
|
protected |
Definition at line 188 of file kEpsilonPhitF.H.
|
protected |
Definition at line 189 of file kEpsilonPhitF.H.
|
protected |
Definition at line 190 of file kEpsilonPhitF.H.