Go to the documentation of this file.
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();
dimensionedScalar sigmaK_
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
volScalarField epsilon_
Turbulent kinetic energy dissipation rate [m2/s3].
virtual void correct()
Solve the transport equations and correct the turbulent viscosity.
A class for handling words, derived from Foam::string.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon (LUU:Eq. 4)
A class for managing temporary objects.
BasicTurbulenceModel::alphaField alphaField
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField f_
Elliptic relaxation factor [1/s].
static const word propertiesName
Default name of the turbulence properties dictionary.
dimensionedScalar sigmaEps_
BasicTurbulenceModel::rhoField rhoField
volScalarField T_
Turbulent time scale [s].
tmp< volScalarField > Ls() const
Return the turbulent length scale, L.
BasicTurbulenceModel::transportModel transportModel
virtual void correctNut()
Update nut with the latest available k, phit, and T.
virtual tmp< volScalarField > phit() const
Return the normalised wall-normal fluctuating velocity scale field.
TypeName("kEpsilonPhitF")
Runtime type information.
dimensionedScalar Ceps1b_
dimensionedScalar phitMin_
virtual tmp< volScalarField > epsilon() const
Return the turbulent kinetic energy dissipation rate field.
dimensionedScalar sigmaPhit_
dimensionedScalar Ceps1a_
tmp< volScalarField > DphitEff() const
Return the effective diffusivity for phit (LUU:Eq. 17)
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k (LUU:Eq. 3)
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy field.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Eddy viscosity turbulence model base class.
virtual ~kEpsilonPhitF()=default
Destructor.
volScalarField k_
Turbulent kinetic energy [m2/s2].
volScalarField phit_
Normalised wall-normal fluctuating velocity scale [-].
dimensionedScalar Ceps1c_
The k-epsilon-phit-f turbulence closure model for incompressible and compressible flows.
virtual bool read()
Re-read model coefficients if they have changed.
virtual tmp< volScalarField > f() const
Return the elliptic relaxation factor field.
tmp< volScalarField > Ts() const
Return the turbulent time scale, T.