Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
35namespace Foam
36{
37namespace regionModels
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
46
47const Enum
48<
50>
52{
53 { frictionMethodType::mquadraticProfile, "quadraticProfile" },
54 { frictionMethodType::mlinearProfile, "linearProfile" },
55 { frictionMethodType::mDarcyWeisbach, "DarcyWeisbach" },
56 { frictionMethodType::mManningStrickler, "ManningStrickler" }
57};
58
59
60const Enum
61<
63>
65{
67 { shearMethodType::mwallFunction, "wallFunction" }
68};
69
70
71// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
72
73filmTurbulenceModel::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 CwPtr_(nullptr),
87 dwfPtr_(nullptr)
88{
89 if (rhoName_ == "rhoInf")
90 {
91 rhoRef_ = dict_.get<scalar>("rhoInf");
92 }
93
94 const auto& regionMesh = film_.regionMesh();
95 CwPtr_.reset
96 (
98 (
100 (
102 regionMesh.time().timeName(),
103 regionMesh.thisDb(),
107 ),
108 regionMesh,
110 )
111 );
112
113 auto* dwfNamePtr = dict_.findEntry("DarcyWeisbachField");
114 if (dwfNamePtr)
115 {
116 dwfPtr_.reset
117 (
119 (
121 (
122 word(dwfNamePtr->stream()),
123 regionMesh.time().timeName(),
124 regionMesh.thisDb(),
127 ),
128 regionMesh
129 )
130 );
131 }
132}
133
134
135// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
138{
139 return film_;
140}
141
142
144{
145 auto& Cw = CwPtr_.ref();
146
147 switch (method_)
148 {
150 {
151 const scalarField& mu = film_.mu().primitiveField();
152 const scalarField& h = film_.h().primitiveField();
154
155 const scalar h0 = film_.h0().value();
156
157 Cw.primitiveFieldRef() = 3*mu/((h + h0)*rho);
158 Cw.clamp_max(5000.0);
159
160 break;
161 }
162 case mlinearProfile:
163 {
164 const scalarField& mu = film_.mu().primitiveField();
165 const scalarField& h = film_.h().primitiveField();
167
168 const scalar h0 = film_.h0().value();
169
170 Cw.primitiveFieldRef() = 2*mu/((h + h0)*rho);
171
172 break;
173 }
174 case mDarcyWeisbach:
175 {
176 const auto& primaryMesh = film_.primaryMesh();
177 const auto& g = meshObjects::gravity::New(primaryMesh.time());
178 const auto magG(mag(g.value()));
179
180 const vectorField& Uf = film_.Uf().primitiveField();
182
183 auto& Cwp = Cw.primitiveFieldRef();
184
185 auto* dwfNamePtr = dict_.findEntry("DarcyWeisbachField");
186 if (dwfNamePtr)
187 {
188 Cwp = dwfPtr_.ref().primitiveField()*magG*mag(Uf)/rho;
189 }
190 else
191 {
192 const scalar Cf = dict_.get<scalar>("DarcyWeisbach");
193 Cwp = Cf*magG*mag(Uf)/rho;
194 }
195
196 break;
197 }
199 {
200 const auto& primaryMesh = film_.primaryMesh();
201 const auto& g = meshObjects::gravity::New(primaryMesh.time());
202
203 const vectorField& Uf = film_.Uf().primitiveField();
204 const scalarField& h = film_.h().primitiveField();
205 const scalar h0 = film_.h0().value();
206
207 const scalar n = dict_.get<scalar>("n");
208
209 Cw.primitiveFieldRef() =
210 sqr(n)*mag(g.value())*mag(Uf)/cbrt(h + h0);
211
212 break;
213 }
214 default:
215 {
217 << "Unimplemented method "
219 << "Please set 'frictionMethod' to one of "
220 << flatOutput(frictionMethodTypeNames_.sortedToc()) << nl
221 << exit(FatalError);
222
223 break;
225 }
226
227 return CwPtr_();
228}
229
230
232(
234) const
235{
236 auto tshearStress =
237 tmp<faVectorMatrix>::New(U, sqr(U.dimensions())*sqr(dimLength));
238
239 switch (shearMethod_)
240 {
241 case msimple:
242 {
244
245 const dimensionedScalar Cf
246 (
247 "Cf",
249 dict_.get<scalar>("Cf")
250 );
251
252 tshearStress.ref() += - fam::Sp(Cf, U) + Cf*Up();
253
254 break;
255 }
256 case mwallFunction:
257 {
258 tmp<volSymmTensorField> tdevRhoReff = devRhoReff();
259
262
263 const volSymmTensorField::Boundary& devRhoReffb
264 = tdevRhoReff().boundaryField();
265
266 // The polyPatch/local-face for each of the faceLabels
267 const List<labelPair>& patchFaces
269
270 // All referenced polyPatches (== primaryPatchIDs())
271 const labelList& patchIds
273
274 // Values per patch
275 PtrMap<vectorField> patchFields(2*patchIds.size());
276
277 for (const label patchi : patchIds)
278 {
279 patchFields.set
280 (
281 patchi,
282 Sfb[patchi] & devRhoReffb[patchi]
283 );
284 }
285
286 vectorField afT(patchFaces.size(), Zero);
287
288 const vectorField& nHat =
289 film_.regionMesh().faceAreaNormals().internalField();
290
291 forAll(patchFaces, i)
292 {
293 const label patchi = patchFaces[i].first();
294 const label facei = patchFaces[i].second();
295
296 const auto* pfld = patchFields.get(patchi);
297
298 if (pfld)
299 {
300 afT[i] = (*pfld)[facei];
301 // Subtract normal component
302 afT[i].removeCollinear(nHat[i]);
303 }
304 }
305
306 auto taForce = areaVectorField::New
307 (
308 "taForce",
310 film_.regionMesh(),
312 );
313 vectorField& aForce = taForce.ref().primitiveFieldRef();
314
315 const DimensionedField<scalar, areaMesh>& magSf =
316 film_.regionMesh().S();
317
318 aForce = afT/(film_.rho().primitiveField()*magSf);
319
320 tshearStress.ref() += taForce();
321
322 if (film_.regionMesh().time().writeTime())
323 {
324 taForce().write();
325 }
326
327 break;
329 }
330
331 return tshearStress;
332}
333
334
336{
337 typedef compressible::turbulenceModel cmpTurbModel;
338 typedef incompressible::turbulenceModel icoTurbModel;
339
340 const fvMesh& m = film_.primaryMesh();
341
342 const auto& U = m.lookupObject<volVectorField>(film_.UName());
343
344 if (m.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
345 {
346 const auto& turb =
347 m.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
348
349 return turb.devRhoReff();
350 }
351 else if (m.foundObject<icoTurbModel>(icoTurbModel::propertiesName))
352 {
353 const auto& turb =
354 m.lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
355
356 return rho()*turb.devReff();
357 }
358 else if (m.foundObject<fluidThermo>(fluidThermo::dictName))
359 {
360 const auto& thermo =
361 m.lookupObject<fluidThermo>(fluidThermo::dictName);
362
363 return -thermo.mu()*devTwoSymm(fvc::grad(U));
364 }
365 else if (m.foundObject<transportModel>("transportProperties"))
366 {
367 const auto& laminarT =
368 m.lookupObject<transportModel>("transportProperties");
369
370 return -rho()*laminarT.nu()*devTwoSymm(fvc::grad(U));
371 }
372 else if (m.foundObject<dictionary>("transportProperties"))
373 {
374 const auto& transportProperties =
375 m.lookupObject<dictionary>("transportProperties");
376
378
379 return -rho()*nu*devTwoSymm(fvc::grad(U));
380 }
381 else
382 {
384 << "No valid model for viscous stress calculation"
386
387 return nullptr;
388 }
389}
390
391
393{
394 const fvMesh& mesh = film_.primaryMesh();
395
396 if (rhoName_ == "rhoInf")
397 {
399 (
400 "rho",
402 mesh,
404 );
405 }
406
407 return mesh.lookupObject<volScalarField>(rhoName_);
408}
409
410
411// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
412
413} // End namespace areaSurfaceFilmModels
414} // End namespace regionModels
415} // End namespace Foam
416
417// ************************************************************************* //
label n
const uniformDimensionedVectorField & g
compressible::turbulenceModel & turb
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
static tmp< GeometricField< vector, faPatchField, areaMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=faPatchField< vector >::calculatedType())
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const T * get(const Key &key) const
Return const pointer associated with given entry or a nullptr if the key does not exist in the table.
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A HashTable of pointers to objects of type <T> with a label key.
Definition PtrMap.H:49
T & first()
Access first element of the list, position [0].
Definition UList.H:957
bool get(const label i) const
Return bool value at specified position, always false for out-of-range access.
Definition UList.H:868
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static const word dictName
The dictionary name ("thermophysicalProperties").
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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 entry * findEntry(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition dictionaryI.H:84
const Type & value() const noexcept
Return const reference to value.
const List< labelPair > & whichPatchFaces() const
The polyPatch/local-face for each faceLabels().
Definition faMeshI.H:138
const labelList & whichPolyPatches() const
The polyPatches related to the areaMesh, in sorted order.
Definition faMeshI.H:128
Fundamental fluid thermodynamic properties.
Definition fluidThermo.H:52
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const surfaceVectorField & Sf() const
Return cell face area vectors.
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
bool foundObject(const word &name, const bool recursive=false) const
Contains the named Type?
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...
virtual tmp< volScalarField > mu() const
Dynamic viscosity of mixture [kg/m/s].
Definition psiThermo.C:175
scalar rhoRef_
Reference density needed for incompressible calculations.
autoPtr< areaScalarField > dwfPtr_
Darcy-Weisbach model field.
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.
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
const liquidFilmBase & film_
Reference to liquidFilmBase.
virtual tmp< areaScalarField > Cw() const
Return the wall film surface friction.
const dimensionedScalar & h0() const noexcept
Return h0.
tmp< areaVectorField > Up() const
Primary region velocity at film hight. Assume the film to be.
virtual const areaScalarField & mu() const =0
Access const reference mu.
virtual const areaScalarField & rho() const =0
Access const reference rho.
const areaVectorField & Uf() const noexcept
Access const reference Uf.
const areaScalarField & h() const noexcept
Access const reference h.
Base class for liquid-film models.
const faMesh & regionMesh() const
Return the region mesh database.
const fvMesh & primaryMesh() const noexcept
Return the reference to the primary mesh database.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Base-class for all transport models used by the incompressible turbulence models.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
labelList patchIds
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
zeroField Sp(const Foam::zero, const GeometricField< Type, faPatchField, areaMesh > &)
A no-op source.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Namespace for OpenFOAM.
const dimensionSet dimViscosity
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
const dimensionSet dimVelocity
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
GeometricField< vector, faPatchField, areaMesh > areaVectorField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
scalar h0
volScalarField & nu
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & h
psiReactionThermo & thermo
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299