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 
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 = areaScalarField::New
105  (
106  "tCw",
108  film_.regionMesh(),
110  );
111  auto& Cw = tCw.ref();
112 
113 
114  switch (method_)
115  {
116  case mquadraticProfile:
117  {
118  const scalarField& mu = film_.mu().primitiveField();
119  const scalarField& h = film_.h().primitiveField();
120  const scalarField& rho = film_.rho().primitiveField();
121 
122  const scalar h0 = film_.h0().value();
123 
124  Cw.primitiveFieldRef() = 3*mu/((h + h0)*rho);
125  Cw.clamp_max(5000.0);
126 
127  break;
128  }
129  case mlinearProfile:
130  {
131  const scalarField& mu = film_.mu().primitiveField();
132  const scalarField& h = film_.h().primitiveField();
133  const scalarField& rho = film_.rho().primitiveField();
134 
135  const scalar h0 = film_.h0().value();
136 
137  Cw.primitiveFieldRef() = 2*mu/((h + h0)*rho);
138 
139  break;
140  }
141  case mDarcyWeisbach:
142  {
145  const vectorField& Uf = film_.Uf().primitiveField();
146  const scalarField& rho = film_.rho().primitiveField();
147 
148  const scalar Cf = dict_.get<scalar>("DarcyWeisbach");
149 
150  Cw.primitiveFieldRef() = Cf*mag(g.value())*mag(Uf)/rho;
151 
152  break;
153  }
154  case mManningStrickler:
155  {
158 
159  const vectorField& Uf = film_.Uf().primitiveField();
160  const scalarField& h = film_.h().primitiveField();
161  const scalar h0 = film_.h0().value();
162 
163  const scalar n = dict_.get<scalar>("n");
164 
165  Cw.primitiveFieldRef() =
166  sqr(n)*mag(g.value())*mag(Uf)/cbrt(h + h0);
167 
168  break;
169  }
170  default:
171  {
173  << "Unimplemented method "
175  << "Please set 'frictionMethod' to one of "
177  << exit(FatalError);
178 
179  break;
180  }
181  }
182 
183  return tCw;
184 }
185 
186 
188 (
190 ) const
191 {
192  auto tshearStress =
193  tmp<faVectorMatrix>::New(U, sqr(U.dimensions())*sqr(dimLength));
194 
195  switch (shearMethod_)
196  {
197  case msimple:
198  {
200 
201  const dimensionedScalar Cf
202  (
203  "Cf",
204  dimVelocity,
205  dict_.get<scalar>("Cf")
206  );
207 
208  tshearStress.ref() += - fam::Sp(Cf, U) + Cf*Up();
209 
210  break;
211  }
212  case mwallFunction:
213  {
214  tmp<volSymmTensorField> tdevRhoReff = devRhoReff();
215 
218 
219  const volSymmTensorField::Boundary& devRhoReffb
220  = tdevRhoReff().boundaryField();
221 
222  // The polyPatch/local-face for each of the faceLabels
223  const List<labelPair>& patchFaces
225 
226  // All referenced polyPatches (== primaryPatchIDs())
227  const labelList& patchIds
229 
230  // Values per patch
231  PtrMap<vectorField> patchFields(2*patchIds.size());
232 
233  for (const label patchi : patchIds)
234  {
235  patchFields.set
236  (
237  patchi,
238  Sfb[patchi] & devRhoReffb[patchi]
239  );
240  }
241 
242  vectorField afT(patchFaces.size(), Zero);
243 
244  const vectorField& nHat =
246 
247  forAll(patchFaces, i)
248  {
249  const label patchi = patchFaces[i].first();
250  const label facei = patchFaces[i].second();
251 
252  const auto* pfld = patchFields.get(patchi);
253 
254  if (pfld)
255  {
256  afT[i] = (*pfld)[facei];
257  // Subtract normal component
258  afT[i].removeCollinear(nHat[i]);
259  }
260  }
261 
262  auto taForce = areaVectorField::New
263  (
264  "taForce",
266  film_.regionMesh(),
268  );
269  vectorField& aForce = taForce.ref().primitiveFieldRef();
270 
272  film_.regionMesh().S();
273 
274  aForce = afT/(film_.rho().primitiveField()*magSf);
275 
276  tshearStress.ref() += taForce();
277 
278  if (film_.regionMesh().time().writeTime())
279  {
280  taForce().write();
281  }
282 
283  break;
284  }
285  }
286 
287  return tshearStress;
288 }
289 
290 
292 {
293  typedef compressible::turbulenceModel cmpTurbModel;
294  typedef incompressible::turbulenceModel icoTurbModel;
295 
296  const fvMesh& m = film_.primaryMesh();
297 
298  const auto& U = m.lookupObject<volVectorField>(film_.UName());
299 
300  if (m.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
301  {
302  const auto& turb =
303  m.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
304 
305  return turb.devRhoReff();
306  }
307  else if (m.foundObject<icoTurbModel>(icoTurbModel::propertiesName))
308  {
309  const auto& turb =
310  m.lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
311 
312  return rho()*turb.devReff();
313  }
315  {
316  const auto& thermo =
318 
319  return -thermo.mu()*devTwoSymm(fvc::grad(U));
320  }
321  else if (m.foundObject<transportModel>("transportProperties"))
322  {
323  const auto& laminarT =
324  m.lookupObject<transportModel>("transportProperties");
325 
326  return -rho()*laminarT.nu()*devTwoSymm(fvc::grad(U));
327  }
328  else if (m.foundObject<dictionary>("transportProperties"))
329  {
330  const auto& transportProperties =
331  m.lookupObject<dictionary>("transportProperties");
332 
334 
335  return -rho()*nu*devTwoSymm(fvc::grad(U));
336  }
337  else
338  {
340  << "No valid model for viscous stress calculation"
342 
343  return nullptr;
344  }
345 }
346 
347 
349 {
350  const fvMesh& mesh = film_.primaryMesh();
351 
352  if (rhoName_ == "rhoInf")
353  {
354  return volScalarField::New
355  (
356  "rho",
358  mesh,
360  );
361  }
362 
364 }
365 
366 
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
368 
369 } // End namespace areaSurfaceFilmModels
370 } // End namespace regionModels
371 } // End namespace Foam
372 
373 // ************************************************************************* //
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)
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:608
tmp< faVectorMatrix > primaryRegionFriction(areaVectorField &U) const
Return primary region friction.
const Time & time() const
Return reference to time.
Definition: faMesh.C:1022
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
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:862
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:1133
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:791
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:128
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)
Base class for liquid-film models.
#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.
dynamicFvMesh & mesh
defineRunTimeSelectionTable(liquidFilmBase, dictionary)
A class for handling words, derived from Foam::string.
Definition: word.H:63
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 areaScalarField & h() const noexcept
Access const reference h.
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...
static tmp< GeometricField< scalar, faPatchField, areaMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=faPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
const uniformDimensionedVectorField & g
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.
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...
const areaVectorField & Uf() const noexcept
Access const reference Uf.
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:138
A class for managing temporary objects.
Definition: HashPtrTable.H:50
scalar h0
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition: faMesh.C:1183
volScalarField & nu
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Do not request registration (bool: false)
Namespace for OpenFOAM.
const dimensionedScalar & h0() const noexcept
Return h0.
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