98 #ifndef kEpsilonPhitF_H 99 #define kEpsilonPhitF_H 115 template<
class BasicTurbulenceModel>
118 public eddyViscosity<RASModel<BasicTurbulenceModel>>
123 kEpsilonPhitF(
const kEpsilonPhitF&) =
delete;
126 void operator=(
const kEpsilonPhitF&) =
delete;
192 typedef typename BasicTurbulenceModel::alphaField
alphaField;
193 typedef typename BasicTurbulenceModel::rhoField
rhoField;
194 typedef typename BasicTurbulenceModel::transportModel
transportModel;
234 this->
nut_/sigmaK_ + this->
nu()
247 this->
nut_/sigmaEps_ + this->
nu()
260 tfld.ref() += this->
nu();
267 virtual tmp<volScalarField>
k()
const virtual void correct()
Solve the transport equations and correct the turbulent viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon (LUU:Eq. 4)
volScalarField phit_
Normalised wall-normal fluctuating velocity scale [-].
volScalarField f_
Elliptic relaxation factor [1/s].
dimensionedScalar Ceps1c_
BasicTurbulenceModel::alphaField alphaField
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
dimensionedScalar sigmaEps_
tmp< volScalarField > Ls() const
Return the turbulent length scale, L.
volScalarField T_
Turbulent time scale [s].
virtual tmp< volScalarField > phit() const
Return the normalised wall-normal fluctuating velocity scale field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
BasicTurbulenceModel::rhoField rhoField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
static const word propertiesName
Default name of the turbulence properties dictionary.
dimensionedScalar phitMin_
A class for handling words, derived from Foam::string.
virtual tmp< volScalarField > epsilon() const
Return the turbulent kinetic energy dissipation rate field.
dimensionedScalar Ceps1a_
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
virtual void correctNut()
Update nut with the latest available k, phit, and T.
BasicTurbulenceModel::transportModel transportModel
TypeName("kEpsilonPhitF")
Runtime type information.
dimensionedScalar Ceps1b_
dimensionedScalar sigmaPhit_
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DphitEff() const
Return the effective diffusivity for phit (LUU:Eq. 17)
tmp< volScalarField > Ts() const
Return the turbulent time scale, T.
The k-epsilon-phit-f turbulence closure model for incompressible and compressible flows...
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k (LUU:Eq. 3)
A class for managing temporary objects.
volScalarField epsilon_
Turbulent kinetic energy dissipation rate [m2/s3].
volScalarField k_
Turbulent kinetic energy [m2/s2].
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > f() const
Return the elliptic relaxation factor field.
dimensionedScalar sigmaK_
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy field.
virtual ~kEpsilonPhitF()=default
Destructor.