Go to the documentation of this file.
37 namespace regionModels
39 namespace areaSurfaceFilmModels
53 { frictionMethodType::mquadraticProfile,
"quadraticProfile" },
54 { frictionMethodType::mlinearProfile,
"linearProfile" },
55 { frictionMethodType::mDarcyWeisbach,
"DarcyWeisbach" },
56 { frictionMethodType::mManningStrickler,
"ManningStrickler" }
66 { shearMethodType::msimple,
"simple" },
67 { shearMethodType::mwallFunction,
"wallFunction" }
73 filmTurbulenceModel::filmTurbulenceModel
75 const word& modelType,
81 dict_(
dict.subDict(modelType +
"Coeffs")),
82 method_(frictionMethodTypeNames_.
get(
"friction", dict_)),
83 shearMethod_(shearMethodTypeNames_.
get(
"shearStress", dict_)),
84 rhoName_(dict_.getOrDefault<
word>(
"rho",
"rho")),
116 auto&
Cw = tCw.ref();
153 const scalar Cf =
dict_.
get<scalar>(
"DarcyWeisbach");
170 Cw.primitiveFieldRef() =
178 <<
"Unimplemented method "
180 <<
"Please set 'frictionMethod' to one of "
224 = tdevRhoReff().boundaryField();
237 fT -= nHat*(fT & nHat);
250 vectorField& aForce = taForce.ref().primitiveFieldRef();
260 tshearStress.
ref() += taForce();
284 if (m.
foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
287 m.
lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
289 return turb.devRhoReff();
291 else if (m.
foundObject<icoTurbModel>(icoTurbModel::propertiesName))
294 m.
lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
307 const auto& laminarT =
324 <<
"No valid model for viscous stress calculation"
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from Foam::string.
faMatrix< vector > faVectorMatrix
const dimensionedScalar mu
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const fvMesh & primaryMesh() const
A class for managing temporary objects.
static constexpr const zero Zero
const dimensionSet dimVelocity
const dimensionSet dimDensity
static const word dictName
const Internal::FieldType & primitiveField() const
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static word timeName(const scalar t, const int precision=precision_)
const shearMethodType shearMethod_
autoPtr< surfaceVectorField > Uf
tmp< Field< Type > > mapToSurface(const typename GeometricField< Type, fvPatchField, volMesh >::Boundary &df) const
Fundamental fluid thermodynamic properties.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
tmp< volScalarField > rho() const
tmp< areaVectorField > Up() const
const Type & value() const
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
const Internal & internalField() const
EnumType get(const word &enumName) const
const volSurfaceMapping & vsm() const
const frictionMethodType method_
bool foundObject(const word &name, const bool recursive=false) const
const areaVectorField & Uf() const
static const Enum< frictionMethodType > frictionMethodTypeNames_
tmp< faMatrix< Type > > Sp(const areaScalarField &sp, const GeometricField< Type, faPatchField, areaMesh > &vf)
const liquidFilmBase & film_
Generic templated field type.
const DimensionedField< scalar, areaMesh > & S() const
bool writeTime() const noexcept
const dimensionedScalar h
const areaVectorField & faceAreaNormals() const
const Time & time() const
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
defineRunTimeSelectionTable(liquidFilmBase, dictionary)
virtual tmp< areaScalarField > Cw() const
const Type & lookupObject(const word &name, const bool recursive=false) const
const dimensionSet dimViscosity
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const faMesh & regionMesh() const
Generic dimensioned Type class.
tmp< faVectorMatrix > primaryRegionFriction(areaVectorField &U) const
Mesh data needed to do the Finite Volume discretisation.
const uniformDimensionedVectorField & g
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Base-class for all transport models used by the incompressible turbulence models.
tmp< volSymmTensorField > devRhoReff() const
virtual const areaScalarField & rho() const =0
#define FatalErrorInFunction
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static const Enum< shearMethodType > shearMethodTypeNames_
defineTypeNameAndDebug(kinematicThinFilm, 0)
const liquidFilmBase & film() const
const dimensionedScalar & h0() const
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static tmp< T > New(Args &&... args)
const areaScalarField & h() const
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
const Time & time() const
List< word > sortedToc() const
Templated abstract base class for single-phase incompressible turbulence models.
static const GeometricField< Type, PatchField, GeoMesh > & null()
dimensionedScalar cbrt(const dimensionedScalar &ds)
Generic GeometricField class.
compressible::turbulenceModel & turb
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const gravity & New(const Time &runTime)
const Boundary & boundaryField() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual const areaScalarField & mu() const =0
Base class for thermal 2D shells.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
const surfaceVectorField & Sf() const