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-------------------------------------------------------------------------------
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
38{
39namespace areaSurfaceFilmModels
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{
66 { shearMethodType::msimple, "simple" },
67 { shearMethodType::mwallFunction, "wallFunction" }
68};
69
70
71// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
72
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 (
108 (
109 "tCw",
112 ),
115 );
116 auto& Cw = tCw.ref();
117
118
119 switch (method_)
120 {
122 {
123 const scalarField& mu = film_.mu().primitiveField();
124 const scalarField& h = film_.h().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();
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();
152
153 const scalar Cf = dict_.get<scalar>("DarcyWeisbach");
154
155 Cw.primitiveFieldRef() = Cf*mag(g.value())*mag(Uf)/rho;
156
157 break;
158 }
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",
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
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 (
242 (
243 "taForce",
246 ),
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"
325 << exit(FatalError);
326
328 }
329}
330
331
333{
334 const fvMesh& m = film_.primaryMesh();
335 if (rhoName_ == "rhoInf")
336 {
338 (
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// ************************************************************************* //
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:61
List< word > sortedToc() const
The sorted list of enum names.
Definition: EnumI.H:70
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const GeometricField< symmTensor, fvPatchField, volMesh > & null()
Return a null geometric field.
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Templated abstract base class for single-phase incompressible turbulence models.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
bool writeTime() const noexcept
True if this is a write time.
Definition: TimeStateI.H:67
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
static const word dictName
Definition: basicThermo.H:256
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
const Type & value() const
Return const reference to value.
const Time & time() const
Return reference to time.
Definition: faMesh.C:673
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Definition: faMesh.C:780
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition: faMesh.C:830
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const surfaceVectorField & Sf() const
Return cell face area vectors.
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
const Type & lookupObject(const word &name, const bool recursive=false) const
scalar rhoRef_
Reference density needed for incompressible calculations.
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)
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 areaVectorField & Uf() const
Access const reference Uf.
const areaScalarField & h() const
Access const reference h.
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 dimensionedScalar & h0() const
Return h0.
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
label patchID() const
Return patch ID.
const volSurfaceMapping & vsm() const
Return mapping between surface and volume fields.
const faMesh & regionMesh() const
Return the region 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:65
T & ref() const
Definition: tmpI.H:227
Base-class for all transport models used by the incompressible turbulence models.
tmp< Field< Type > > mapToSurface(const GeometricBoundaryField< Type, fvPatchField, volMesh > &df) const
Map volume boundary field to surface.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
const volScalarField & mu
autoPtr< surfaceVectorField > Uf
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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:54
Namespace for OpenFOAM.
const dimensionSet dimViscosity
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimVelocity
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:215
faMatrix< vector > faVectorMatrix
Definition: faMatrices.H:53
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
dimensionedScalar cbrt(const dimensionedScalar &ds)
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
scalar h0
volScalarField & nu
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & h
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))