133#ifndef Foam_functionObjects_turbulenceFields_H
134#define Foam_functionObjects_turbulenceFields_H
233 const word& fieldName,
238 template<
class Model>
242 template<
class Model>
246 template<
class Model>
287 virtual bool write();
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Generic GeometricField class.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
Computes various turbulence-related quantities that are not typically output during calculations,...
static const word modelName_
Name of the turbulence properties dictionary.
word prefix_
Name of output-field prefix.
void initialise()
Unset duplicate fields already registered by other function objects.
turbulenceFields(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
void operator=(const turbulenceFields &)=delete
No copy assignment.
bool compressible()
Return true if compressible turbulence model is identified.
bool initialised_
Flag to track initialisation.
tmp< volScalarField > I(const Model &model) const
Return turbulence intensity, I, calculated from k and U.
tmp< volScalarField > nuTilda(const Model &model) const
Return nuTilda calculated from k and omega.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
static const Enum< incompressibleField > incompressibleFieldNames_
Names for incompressibleField turbulence fields.
compressibleField
Options for the turbulence fields (compressible).
@ cfK
"Turbulent kinetic energy"
@ cfAlphat
"Turbulence thermal diffusivity"
@ cfL
"Integral-length/Mixing-length scale"
@ cfNuTilda
"Modified turbulent viscosity"
@ cfDevRhoReff
"Divergence of the Reynolds stress"
@ cffd
"DES model shielding function"
@ cfLESRegion
"DES model LES region indicator field"
@ cfMuEff
"Effective turbulent dynamic viscosity"
@ cfAlphaEff
"Effective turbulence thermal diffusivity"
@ cfI
"Turbulence intensity"
@ cfOmega
"Specific dissipation rate"
@ cfMut
"Turbulent dynamic viscosity"
@ cfR
"Reynolds stress tensor"
@ cfEpsilon
"Turbulent kinetic energy dissipation rate"
incompressibleField
Options for the turbulence fields (incompressible).
@ ifK
"Turbulent kinetic energy"
@ ifI
"Turbulence intensity"
@ ifEpsilon
"Turbulent kinetic energy dissipation rate"
@ ifNut
"Turbulent viscosity"
@ iffd
"DES model shielding function"
@ ifL
"Integral-length/Mixing-length scale"
@ ifLESRegion
"DES model LES region indicator field"
@ ifNuTilda
"Modified turbulent viscosity"
@ ifR
"Reynolds stress tensor"
@ ifOmega
"Specific dissipation rate"
@ ifDevReff
"Deviatoric part of the effective Reynolds stress"
@ ifNuEff
"Effective turbulent viscosity"
void processField(const word &fieldName, const tmp< GeometricField< Type, fvPatchField, volMesh > > &tvalue)
Process the turbulence field.
turbulenceFields(const turbulenceFields &)=delete
No copy construct.
TypeName("turbulenceFields")
Runtime type information.
wordHashSet fieldSet_
Fields to load.
virtual bool execute()
Execute the function-object operations.
tmp< volScalarField > L(const Model &model) const
Return integral length scale, L, calculated from k and epsilon.
virtual ~turbulenceFields()=default
Destructor.
static const Enum< compressibleField > compressibleFieldNames_
Names for compressibleField turbulence fields.
virtual bool write()
Write the function-object results.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
static const Identity< scalar > I
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
const vector L(dict.get< vector >("L"))
Forwards and collection of common volume field types.