filmTurbulenceModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2020-2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "filmTurbulenceModel.H"
29 #include "gravityMeshObject.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace regionModels
38 {
39 namespace areaSurfaceFilmModels
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(filmTurbulenceModel, 0);
45 defineRunTimeSelectionTable(filmTurbulenceModel, dictionary);
46 
47 const Enum
48 <
50 >
52 {
53  { frictionMethodType::mquadraticProfile, "quadraticProfile" },
54  { frictionMethodType::mlinearProfile, "linearProfile" },
55  { frictionMethodType::mDarcyWeisbach, "DarcyWeisbach" },
56  { frictionMethodType::mManningStrickler, "ManningStrickler" }
57 };
58 
59 
60 const Enum
61 <
63 >
65 {
66  { shearMethodType::msimple, "simple" },
67  { shearMethodType::mwallFunction, "wallFunction" }
68 };
69 
70 
71 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
72 
73 filmTurbulenceModel::filmTurbulenceModel
74 (
75  const word& modelType,
76  liquidFilmBase& film,
77  const dictionary& dict
78 )
79 :
80  film_(film),
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")),
85  rhoRef_(VGREAT)
86 {
87  if (rhoName_ == "rhoInf")
88  {
89  rhoRef_ = dict_.get<scalar>("rhoInf");
90  }
91 }
92 
93 
94 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
95 
97 {
98  return film_;
99 }
100 
101 
103 {
104  auto tCw =
106  (
107  IOobject
108  (
109  "tCw",
112  ),
113  film_.regionMesh(),
115  );
116  auto& Cw = tCw.ref();
117 
118 
119  switch (method_)
120  {
121  case mquadraticProfile:
122  {
123  const scalarField& mu = film_.mu().primitiveField();
124  const scalarField& h = film_.h().primitiveField();
125  const scalarField& rho = film_.rho().primitiveField();
126 
127  const scalar h0 = film_.h0().value();
128 
129  Cw.primitiveFieldRef() = 3*mu/((h + h0)*rho);
130  Cw.min(5000.0);
131 
132  break;
133  }
134  case mlinearProfile:
135  {
136  const scalarField& mu = film_.mu().primitiveField();
137  const scalarField& h = film_.h().primitiveField();
138  const scalarField& rho = film_.rho().primitiveField();
139 
140  const scalar h0 = film_.h0().value();
141 
142  Cw.primitiveFieldRef() = 2*mu/((h + h0)*rho);
143 
144  break;
145  }
146  case mDarcyWeisbach:
147  {
150  const vectorField& Uf = film_.Uf().primitiveField();
151  const scalarField& rho = film_.rho().primitiveField();
152 
153  const scalar Cf = dict_.get<scalar>("DarcyWeisbach");
154 
155  Cw.primitiveFieldRef() = Cf*mag(g.value())*mag(Uf)/rho;
156 
157  break;
158  }
159  case mManningStrickler:
160  {
163 
164  const vectorField& Uf = film_.Uf().primitiveField();
165  const scalarField& h = film_.h().primitiveField();
166  const scalar h0 = film_.h0().value();
167 
168  const scalar n = dict_.get<scalar>("n");
169 
170  Cw.primitiveFieldRef() =
171  sqr(n)*mag(g.value())*mag(Uf)/cbrt(h + h0);
172 
173  break;
174  }
175  default:
176  {
178  << "Unimplemented method "
180  << "Please set 'frictionMethod' to one of "
182  << exit(FatalError);
183 
184  break;
185  }
186  }
187 
188  return tCw;
189 }
190 
191 
193 (
195 ) const
196 {
197  tmp<faVectorMatrix> tshearStress
198  (
199  new faVectorMatrix(U, sqr(U.dimensions())*sqr(dimLength))
200  );
201 
202  switch (shearMethod_)
203  {
204  case msimple:
205  {
207 
208  const dimensionedScalar Cf
209  (
210  "Cf",
211  dimVelocity,
212  dict_.get<scalar>("Cf")
213  );
214 
215  tshearStress.ref() += - fam::Sp(Cf, U) + Cf*Up();
216 
217  break;
218  }
219  case mwallFunction:
220  {
221  tmp<volSymmTensorField> tdevRhoReff = devRhoReff();
222 
223  const volSymmTensorField::Boundary& devRhoReffb
224  = tdevRhoReff().boundaryField();
225 
226  const label patchi = film_.patchID();
227 
228  const surfaceVectorField::Boundary& Sfb =
230 
231  vectorField fT(Sfb[patchi] & devRhoReffb[patchi]);
232 
233  const vectorField& nHat =
235 
236  // Substract normal component
237  fT -= nHat*(fT & nHat);
238 
239  auto taForce = tmp<areaVectorField>::New
240  (
241  IOobject
242  (
243  "taForce",
246  ),
247  film_.regionMesh(),
249  );
250  vectorField& aForce = taForce.ref().primitiveFieldRef();
251 
252  // Map ft to surface
253  const vectorField afT(film_.vsm().mapToSurface(fT));
254 
256  film_.regionMesh().S();
257 
258  aForce = afT/(film_.rho().primitiveField()*magSf);
259 
260  tshearStress.ref() += taForce();
261 
262  if (film_.regionMesh().time().writeTime())
263  {
264  taForce().write();
265  }
266 
267  break;
268  }
269  }
270 
271  return tshearStress;
272 }
273 
274 
276 {
277  typedef compressible::turbulenceModel cmpTurbModel;
278  typedef incompressible::turbulenceModel icoTurbModel;
279 
280  const fvMesh& m = film_.primaryMesh();
281 
282  const auto& U = m.lookupObject<volVectorField>(film_.UName());
283 
284  if (m.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
285  {
286  const auto& turb =
287  m.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
288 
289  return turb.devRhoReff();
290  }
291  else if (m.foundObject<icoTurbModel>(icoTurbModel::propertiesName))
292  {
293  const auto& turb =
294  m.lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
295 
296  return rho()*turb.devReff();
297  }
299  {
300  const auto& thermo =
302 
303  return -thermo.mu()*dev(twoSymm(fvc::grad(U)));
304  }
305  else if (m.foundObject<transportModel>("transportProperties"))
306  {
307  const auto& laminarT =
308  m.lookupObject<transportModel>("transportProperties");
309 
310  return -rho()*laminarT.nu()*dev(twoSymm(fvc::grad(U)));
311  }
312  else if (m.foundObject<dictionary>("transportProperties"))
313  {
314  const auto& transportProperties =
315  m.lookupObject<dictionary>("transportProperties");
316 
318 
319  return -rho()*nu*dev(twoSymm(fvc::grad(U)));
320  }
321  else
322  {
324  << "No valid model for viscous stress calculation"
326 
327  return volSymmTensorField::null();
328  }
329 }
330 
331 
333 {
334  const fvMesh& m = film_.primaryMesh();
335  if (rhoName_ == "rhoInf")
336  {
338  (
339  IOobject
340  (
341  "rho",
342  m.time().timeName(),
343  m
344  ),
345  m,
347  );
348  }
349 
351 }
352 
353 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
354 
355 } // End namespace areaSurfaceFilmModels
356 } // End namespace regionModels
357 } // End namespace Foam
358 
359 // ************************************************************************* //
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::mManningStrickler
@ mManningStrickler
Definition: filmTurbulenceModel.H:77
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase
Definition: liquidFilmBase.H:56
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::shearMethodType
shearMethodType
Definition: filmTurbulenceModel.H:81
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:53
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::faVectorMatrix
faMatrix< vector > faVectorMatrix
Definition: faMatrices.H:47
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::rhoName_
word rhoName_
Definition: filmTurbulenceModel.H:111
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Definition: createFieldRefs.H:4
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::regionModels::regionFaModel::primaryMesh
const fvMesh & primaryMesh() const
Definition: regionFaModelI.H:26
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
Foam::basicThermo::dictName
static const word dictName
Definition: basicThermo.H:252
h0
scalar h0
Definition: readInitialConditions.H:78
filmTurbulenceModel.H
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Definition: GeometricFieldI.H:46
transportProperties
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
turbulentTransportModel.H
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::shearMethod_
const shearMethodType shearMethod_
Definition: filmTurbulenceModel.H:108
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:27
gravityMeshObject.H
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::frictionMethodType
frictionMethodType
Definition: filmTurbulenceModel.H:72
Foam::volSurfaceMapping::mapToSurface
tmp< Field< Type > > mapToSurface(const typename GeometricField< Type, fvPatchField, volMesh >::Boundary &df) const
Foam::fluidThermo
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:48
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::rho
tmp< volScalarField > rho() const
Definition: filmTurbulenceModel.C:325
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::Up
tmp< areaVectorField > Up() const
Definition: liquidFilmBase.C:344
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::UName
word UName() const
Definition: liquidFilmBase.C:548
Foam::dimensioned::value
const Type & value() const
Definition: dimensionedType.C:427
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:100
Foam::GeometricField::internalField
const Internal & internalField() const
Definition: GeometricFieldI.H:36
Foam::Enum::get
EnumType get(const word &enumName) const
Definition: Enum.C:68
Foam::regionModels::regionFaModel::vsm
const volSurfaceMapping & vsm() const
Definition: regionFaModel.C:103
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::method_
const frictionMethodType method_
Definition: filmTurbulenceModel.H:105
Foam::objectRegistry::foundObject
bool foundObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:372
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::Uf
const areaVectorField & Uf() const
Definition: liquidFilmBase.C:494
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::rhoRef_
scalar rhoRef_
Definition: filmTurbulenceModel.H:114
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::frictionMethodTypeNames_
static const Enum< frictionMethodType > frictionMethodTypeNames_
Definition: filmTurbulenceModel.H:96
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::fam::Sp
tmp< faMatrix< Type > > Sp(const areaScalarField &sp, const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famSup.C:79
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::film_
const liquidFilmBase & film_
Definition: filmTurbulenceModel.H:93
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:220
Foam::UniformDimensionedField< vector >
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::faMesh::S
const DimensionedField< scalar, areaMesh > & S() const
Definition: faMesh.C:625
Foam::TimeState::writeTime
bool writeTime() const noexcept
Definition: TimeStateI.H:60
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::dict_
const dictionary dict_
Definition: filmTurbulenceModel.H:102
Foam::constant::universal::h
const dimensionedScalar h
Definition: setRegionSolidFields.H:33
Foam::faMesh::faceAreaNormals
const areaVectorField & faceAreaNormals() const
Definition: faMesh.C:675
Foam::faMesh::time
const Time & time() const
Definition: faMesh.C:539
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::mDarcyWeisbach
@ mDarcyWeisbach
Definition: filmTurbulenceModel.H:76
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:453
Foam::regionModels::areaSurfaceFilmModels::defineRunTimeSelectionTable
defineRunTimeSelectionTable(liquidFilmBase, dictionary)
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::Cw
virtual tmp< areaScalarField > Cw() const
Definition: filmTurbulenceModel.C:95
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:427
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::msimple
@ msimple
Definition: filmTurbulenceModel.H:83
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dimViscosity
const dimensionSet dimViscosity
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::mlinearProfile
@ mlinearProfile
Definition: filmTurbulenceModel.H:75
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::mwallFunction
@ mwallFunction
Definition: filmTurbulenceModel.H:84
Foam::regionModels::regionFaModel::regionMesh
const faMesh & regionMesh() const
Definition: regionFaModelI.H:56
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:36
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::primaryRegionFriction
tmp< faVectorMatrix > primaryRegionFriction(areaVectorField &U) const
Definition: filmTurbulenceModel.C:186
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
Foam
Definition: atmBoundaryLayer.C:26
Foam::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Definition: FlatOutput.H:217
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::transportModel
Base-class for all transport models used by the incompressible turbulence models.
Definition: transportModel.H:49
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::devRhoReff
tmp< volSymmTensorField > devRhoReff() const
Definition: filmTurbulenceModel.C:268
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::rho
virtual const areaScalarField & rho() const =0
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:44
Foam::GeometricField::Boundary
Definition: GeometricField.H:111
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::shearMethodTypeNames_
static const Enum< shearMethodType > shearMethodTypeNames_
Definition: filmTurbulenceModel.H:99
Foam::regionModels::areaSurfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicThinFilm, 0)
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::film
const liquidFilmBase & film() const
Definition: filmTurbulenceModel.C:89
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::h0
const dimensionedScalar & h0() const
Definition: liquidFilmBase.C:530
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:41
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::tmp::New
static tmp< T > New(Args &&... args)
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::h
const areaScalarField & h() const
Definition: liquidFilmBase.C:512
Foam::regionModels::regionFaModel::patchID
label patchID() const
Definition: regionFaModelI.H:131
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::fvMesh::time
const Time & time() const
Definition: fvMesh.H:276
Foam::Enum::sortedToc
List< word > sortedToc() const
Definition: EnumI.H:63
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:51
Foam::GeometricField::null
static const GeometricField< Type, PatchField, GeoMesh > & null()
Definition: GeometricFieldI.H:25
Foam::regionModels::areaSurfaceFilmModels::filmTurbulenceModel::mquadraticProfile
@ mquadraticProfile
Definition: filmTurbulenceModel.H:74
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:148
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:88
turbulentFluidThermoModel.H
Foam::meshObjects::gravity::New
static const gravity & New(const Time &runTime)
Definition: gravityMeshObject.H:89
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Definition: GeometricFieldI.H:55
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:50
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase::mu
virtual const areaScalarField & mu() const =0
liquidFilmBase
Base class for thermal 2D shells.
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:99
Foam::fvMesh::Sf
const surfaceVectorField & Sf() const
Definition: fvMeshGeometry.C:312