101#ifndef regionFaModels_filmTurbulenceModel_H
102#define regionFaModels_filmTurbulenceModel_H
125 filmTurbulenceModel(
const filmTurbulenceModel&) =
delete;
128 void operator=(
const filmTurbulenceModel&) =
delete;
157 const liquidFilmBase&
film_;
213 const word& modelType,
Forwards and collection of common area field types.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Base class for film turbulence models.
scalar rhoRef_
Reference density needed for incompressible calculations.
autoPtr< areaScalarField > dwfPtr_
Darcy-Weisbach model field.
virtual tmp< areaScalarField > mut() const =0
Return the film turbulence viscosity.
frictionMethodType
Options for the friction models.
const liquidFilmBase & film() const
Return film.
static const Enum< frictionMethodType > frictionMethodTypeNames_
Names for friction models.
tmp< faVectorMatrix > primaryRegionFriction(areaVectorField &U) const
Return primary region friction.
static const Enum< shearMethodType > shearMethodTypeNames_
Names for shear stress models.
tmp< volSymmTensorField > devRhoReff() const
Return the effective viscous stress (laminar + turbulent).
autoPtr< areaScalarField > CwPtr_
Wall film-surface friction field.
virtual tmp< faVectorMatrix > Su(areaVectorField &U) const =0
Return the source for the film momentum equation.
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
virtual ~filmTurbulenceModel()=default
Destructor.
TypeName("filmTurbulenceModel")
Runtime type information.
const frictionMethodType method_
Friction model.
const liquidFilmBase & film_
Reference to liquidFilmBase.
virtual void correct()=0
Correct/update the model.
const shearMethodType shearMethod_
Shear-stress model.
static autoPtr< filmTurbulenceModel > New(liquidFilmBase &film, const dictionary &dict)
Return a reference to the selected injection model.
word rhoName_
Name of density field.
const dictionary dict_
Model dictionary.
virtual tmp< areaScalarField > Cw() const
Return the wall film surface friction.
shearMethodType
Options for the shear stress models.
declareRunTimeSelectionTable(autoPtr, filmTurbulenceModel, dictionary,(liquidFilmBase &film, const dictionary &dict),(film, dict))
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
GeometricField< vector, faPatchField, areaMesh > areaVectorField
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.