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-2023 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 #include "PtrMap.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace regionModels
39 {
40 namespace areaSurfaceFilmModels
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 defineTypeNameAndDebug(filmTurbulenceModel, 0);
46 defineRunTimeSelectionTable(filmTurbulenceModel, dictionary);
47 
48 const Enum
49 <
51 >
53 {
54  { frictionMethodType::mquadraticProfile, "quadraticProfile" },
55  { frictionMethodType::mlinearProfile, "linearProfile" },
56  { frictionMethodType::mDarcyWeisbach, "DarcyWeisbach" },
57  { frictionMethodType::mManningStrickler, "ManningStrickler" }
58 };
59 
60 
61 const Enum
62 <
64 >
66 {
67  { shearMethodType::msimple, "simple" },
68  { shearMethodType::mwallFunction, "wallFunction" }
69 };
70 
71 
72 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
73 
74 filmTurbulenceModel::filmTurbulenceModel
75 (
76  const word& modelType,
77  liquidFilmBase& film,
78  const dictionary& dict
79 )
80 :
81  film_(film),
82  dict_(dict.subDict(modelType + "Coeffs")),
83  method_(frictionMethodTypeNames_.get("friction", dict_)),
84  shearMethod_(shearMethodTypeNames_.get("shearStress", dict_)),
85  rhoName_(dict_.getOrDefault<word>("rho", "rho")),
86  rhoRef_(VGREAT)
87 {
88  if (rhoName_ == "rhoInf")
89  {
90  rhoRef_ = dict_.get<scalar>("rhoInf");
91  }
92 }
93 
94 
95 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
96 
98 {
99  return film_;
100 }
101 
102 
104 {
105  auto tCw =
107  (
108  IOobject
109  (
110  "tCw",
113  ),
114  film_.regionMesh(),
116  );
117  auto& Cw = tCw.ref();
118 
119 
120  switch (method_)
121  {
122  case mquadraticProfile:
123  {
124  const scalarField& mu = film_.mu().primitiveField();
125  const scalarField& h = film_.h().primitiveField();
126  const scalarField& rho = film_.rho().primitiveField();
127 
128  const scalar h0 = film_.h0().value();
129 
130  Cw.primitiveFieldRef() = 3*mu/((h + h0)*rho);
131  Cw.clamp_max(5000.0);
132 
133  break;
134  }
135  case mlinearProfile:
136  {
137  const scalarField& mu = film_.mu().primitiveField();
138  const scalarField& h = film_.h().primitiveField();
139  const scalarField& rho = film_.rho().primitiveField();
140 
141  const scalar h0 = film_.h0().value();
142 
143  Cw.primitiveFieldRef() = 2*mu/((h + h0)*rho);
144 
145  break;
146  }
147  case mDarcyWeisbach:
148  {
151  const vectorField& Uf = film_.Uf().primitiveField();
152  const scalarField& rho = film_.rho().primitiveField();
153 
154  const scalar Cf = dict_.get<scalar>("DarcyWeisbach");
155 
156  Cw.primitiveFieldRef() = Cf*mag(g.value())*mag(Uf)/rho;
157 
158  break;
159  }
160  case mManningStrickler:
161  {
164 
165  const vectorField& Uf = film_.Uf().primitiveField();
166  const scalarField& h = film_.h().primitiveField();
167  const scalar h0 = film_.h0().value();
168 
169  const scalar n = dict_.get<scalar>("n");
170 
171  Cw.primitiveFieldRef() =
172  sqr(n)*mag(g.value())*mag(Uf)/cbrt(h + h0);
173 
174  break;
175  }
176  default:
177  {
179  << "Unimplemented method "
181  << "Please set 'frictionMethod' to one of "
183  << exit(FatalError);
184 
185  break;
186  }
187  }
188 
189  return tCw;
190 }
191 
192 
194 (
196 ) const
197 {
198  tmp<faVectorMatrix> tshearStress
199  (
200  new faVectorMatrix(U, sqr(U.dimensions())*sqr(dimLength))
201  );
202 
203  switch (shearMethod_)
204  {
205  case msimple:
206  {
208 
209  const dimensionedScalar Cf
210  (
211  "Cf",
212  dimVelocity,
213  dict_.get<scalar>("Cf")
214  );
215 
216  tshearStress.ref() += - fam::Sp(Cf, U) + Cf*Up();
217 
218  break;
219  }
220  case mwallFunction:
221  {
222  tmp<volSymmTensorField> tdevRhoReff = devRhoReff();
223 
226 
227  const volSymmTensorField::Boundary& devRhoReffb
228  = tdevRhoReff().boundaryField();
229 
230  // The polyPatch/local-face for each of the faceLabels
231  const List<labelPair>& patchFaces
233 
234  // All referenced polyPatches (== primaryPatchIDs())
235  const labelList& patchIds
237 
238  // Values per patch
239  PtrMap<vectorField> patchFields(2*patchIds.size());
240 
241  for (const label patchi : patchIds)
242  {
243  patchFields.set
244  (
245  patchi,
246  Sfb[patchi] & devRhoReffb[patchi]
247  );
248  }
249 
250  vectorField afT(patchFaces.size(), Zero);
251 
252  const vectorField& nHat =
254 
255  forAll(patchFaces, i)
256  {
257  const label patchi = patchFaces[i].first();
258  const label facei = patchFaces[i].second();
259 
260  const auto* pfld = patchFields.get(patchi);
261 
262  if (pfld)
263  {
264  afT[i] = (*pfld)[facei];
265  // Subtract normal component
266  afT[i].removeCollinear(nHat[i]);
267  }
268  }
269 
270  auto taForce = tmp<areaVectorField>::New
271  (
272  IOobject
273  (
274  "taForce",
277  ),
278  film_.regionMesh(),
280  );
281  vectorField& aForce = taForce.ref().primitiveFieldRef();
282 
284  film_.regionMesh().S();
285 
286  aForce = afT/(film_.rho().primitiveField()*magSf);
287 
288  tshearStress.ref() += taForce();
289 
290  if (film_.regionMesh().time().writeTime())
291  {
292  taForce().write();
293  }
294 
295  break;
296  }
297  }
298 
299  return tshearStress;
300 }
301 
302 
304 {
305  typedef compressible::turbulenceModel cmpTurbModel;
306  typedef incompressible::turbulenceModel icoTurbModel;
307 
308  const fvMesh& m = film_.primaryMesh();
309 
310  const auto& U = m.lookupObject<volVectorField>(film_.UName());
311 
312  if (m.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
313  {
314  const auto& turb =
315  m.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
316 
317  return turb.devRhoReff();
318  }
319  else if (m.foundObject<icoTurbModel>(icoTurbModel::propertiesName))
320  {
321  const auto& turb =
322  m.lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
323 
324  return rho()*turb.devReff();
325  }
327  {
328  const auto& thermo =
330 
331  return -thermo.mu()*devTwoSymm(fvc::grad(U));
332  }
333  else if (m.foundObject<transportModel>("transportProperties"))
334  {
335  const auto& laminarT =
336  m.lookupObject<transportModel>("transportProperties");
337 
338  return -rho()*laminarT.nu()*devTwoSymm(fvc::grad(U));
339  }
340  else if (m.foundObject<dictionary>("transportProperties"))
341  {
342  const auto& transportProperties =
343  m.lookupObject<dictionary>("transportProperties");
344 
346 
347  return -rho()*nu*devTwoSymm(fvc::grad(U));
348  }
349  else
350  {
352  << "No valid model for viscous stress calculation"
354 
355  return volSymmTensorField::null();
356  }
357 }
358 
359 
361 {
362  const fvMesh& m = film_.primaryMesh();
363  if (rhoName_ == "rhoInf")
364  {
366  (
367  IOobject
368  (
369  "rho",
370  m.time().timeName(),
371  m
372  ),
373  m,
375  );
376  }
377 
379 }
380 
381 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
382 
383 } // End namespace areaSurfaceFilmModels
384 } // End namespace regionModels
385 } // End namespace Foam
386 
387 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
labelList patchIds
const Type & value() const noexcept
Return const reference to value.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const surfaceVectorField & Sf() const
Return cell face area vectors.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
static const Enum< shearMethodType > shearMethodTypeNames_
Names for shear stress models.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const areaScalarField & h() const
Access const reference h.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tmp< areaVectorField > Up() const
Primary region velocity at film hight. Assume the film to be.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
tmp< faVectorMatrix > primaryRegionFriction(areaVectorField &U) const
Return primary region friction.
const Time & time() const
Return reference to time.
Definition: faMesh.C:791
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static const Enum< frictionMethodType > frictionMethodTypeNames_
Names for friction models.
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
const dimensionSet dimViscosity
compressible::turbulenceModel & turb
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Definition: faMesh.C:914
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type get(const label i) const
Return bool value at specified position, always false for out-of-range access.
Definition: UList.H:782
A HashTable of pointers to objects of type <T> with a label key.
Definition: HashTableFwd.H:42
const labelList & whichPolyPatches() const
The polyPatches related to the areaMesh, in sorted order.
Definition: faMeshI.H:135
faMatrix< vector > faVectorMatrix
Definition: faMatrices.H:47
bool writeTime() const noexcept
True if this is a write interval.
Definition: TimeStateI.H:73
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
defineTypeNameAndDebug(kinematicThinFilm, 0)
virtual const objectRegistry & thisDb() const
Return reference to the mesh database.
Definition: faMesh.C:815
Base class for thermal 2D shells.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Templated abstract base class for single-phase incompressible turbulence models.
tmp< volSymmTensorField > devRhoReff() const
Return the effective viscous stress (laminar + turbulent)
zeroField Sp(const Foam::zero, const GeometricField< Type, faPatchField, areaMesh > &)
A no-op source.
const liquidFilmBase & film_
Reference to liquidFilmBase.
const dimensionedScalar & h0() const
Return h0.
defineRunTimeSelectionTable(liquidFilmBase, dictionary)
A class for handling words, derived from Foam::string.
Definition: word.H:63
const areaVectorField & Uf() const
Access const reference Uf.
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:48
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
dimensionedScalar cbrt(const dimensionedScalar &ds)
const fvMesh & primaryMesh() const noexcept
Return the reference to the primary mesh database.
const faMesh & regionMesh() const
Return the region mesh database.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
autoPtr< surfaceVectorField > Uf
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
const uniformDimensionedVectorField & g
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
const dimensionSet dimDensity
const dimensionedScalar h
Planck constant.
virtual tmp< areaScalarField > Cw() const
Return the wall film surface friction.
static const GeometricField< symmTensor, fvPatchField, volMesh > & null()
Return a null geometric field.
const dimensionedScalar mu
Atomic mass unit.
List< word > sortedToc() const
The sorted list of enum names.
Definition: EnumI.H:63
scalar rhoRef_
Reference density needed for incompressible calculations.
U
Definition: pEqn.H:72
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Base-class for all transport models used by the incompressible turbulence models. ...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual const areaScalarField & rho() const =0
Access const reference rho.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
label n
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
virtual const areaScalarField & mu() const =0
Access const reference mu.
const List< labelPair > & whichPatchFaces() const
The polyPatch/local-face for each faceLabels()
Definition: faMeshI.H:145
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
scalar h0
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition: faMesh.C:964
volScalarField & nu
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Namespace for OpenFOAM.
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity