Loading...
Searching...
No Matches
liquidFilmBase.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-------------------------------------------------------------------------------
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 "liquidFilmBase.H"
29#include "gravityMeshObject.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace regionModels
39{
41{
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
48
49const Foam::word liquidFilmName("liquidFilm");
50
51// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52
54(
55 const word& modelType,
56 const fvMesh& mesh,
57 const dictionary& dict
58)
59:
60 regionFaModel(mesh, liquidFilmName, modelType, dict, true),
62 (
63 this->solution().subDict("PIMPLE").get<bool>("momentumPredictor")
64 ),
66 (
67 this->solution().subDict("PIMPLE").get<label>("nOuterCorr")
68 ),
69 nCorr_(this->solution().subDict("PIMPLE").get<label>("nCorr")),
71 (
72 this->solution().subDict("PIMPLE").get<label>("nFilmCorr")
73 ),
74 h0_("h0", dimLength, 1e-7, dict),
75 deltaWet_("deltaWet", dimLength, 1e-4, dict),
76 UName_(dict.get<word>("U")),
78 pRef_(dict.get<scalar>("pRef")),
79 h_
80 (
81 IOobject
82 (
83 "hf_" + regionName_,
85 regionMesh().thisDb(),
86 IOobject::MUST_READ,
87 IOobject::AUTO_WRITE
88 ),
90 ),
91 Uf_
92 (
93 IOobject
94 (
95 "Uf_" + regionName_,
97 regionMesh().thisDb(),
98 IOobject::MUST_READ,
99 IOobject::AUTO_WRITE
100 ),
101 regionMesh()
102 ),
103 pf_
104 (
105 IOobject
106 (
107 "pf_" + regionName_,
109 regionMesh().thisDb(),
110 IOobject::NO_READ,
111 IOobject::AUTO_WRITE
112 ),
113 regionMesh(),
115 ),
116 ppf_
117 (
118 IOobject
119 (
120 "ppf_" + regionName_,
122 regionMesh().thisDb(),
123 IOobject::NO_READ,
124 IOobject::NO_WRITE
125 ),
126 regionMesh(),
128 ),
129 phif_
130 (
131 IOobject
132 (
133 "phif_" + regionName_,
135 regionMesh().thisDb(),
136 IOobject::READ_IF_PRESENT,
137 IOobject::AUTO_WRITE
138 ),
139 fac::interpolate(Uf_) & regionMesh().Le()
140 ),
141 phi2s_
142 (
143 IOobject
144 (
145 "phi2s_" + regionName_,
147 regionMesh().thisDb(),
148 IOobject::READ_IF_PRESENT,
149 IOobject::AUTO_WRITE
150 ),
151 fac::interpolate(h_*Uf_) & regionMesh().Le()
152 ),
153 gn_
154 (
155 IOobject
156 (
157 "gn",
159 regionMesh().thisDb(),
160 IOobject::NO_READ,
161 IOobject::NO_WRITE
162 ),
163 regionMesh(),
165 ),
166 g_(meshObjects::gravity::New(primaryMesh().time())),
168 (
169 IOobject
170 (
171 "massSource",
173 primaryMesh().thisDb()
174 ),
175 primaryMesh(),
177 ),
179 (
180 IOobject
181 (
182 "momentumSource",
184 primaryMesh().thisDb()
185 ),
186 primaryMesh(),
188 ),
190 (
191 IOobject
192 (
193 "pnSource",
195 primaryMesh().thisDb()
196 ),
197 primaryMesh(),
199 ),
202 (
203 Foam::fa::options::New(primaryMesh(), regionFaModel::areaName())
204 )
205{
207
208 gn_ = g_ & ns;
209
210 if (faOptions_.optionList::empty())
211 {
212 Info<< "No finite area options present for area : "
215}
216
217
218// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
221{}
222
223
224// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
225
227{
228 scalar CoNum = 0.0;
229 scalar velMag = 0.0;
230
231 edgeScalarField SfUfbyDelta
232 (
234 );
235
236 CoNum =
237 max(SfUfbyDelta/regionMesh().magLe()).value()*time().deltaTValue();
238
239 velMag = max(mag(phif_)/regionMesh().magLe()).value();
240
243
244 Info<< "Max film Courant Number: " << CoNum
245 << " Film velocity magnitude: " << velMag << endl;
246
247 return CoNum;
248}
249
250
252{
253 auto tUw = areaVectorField::New
254 (
255 "tUw",
257 regionMesh(),
259 );
260 auto& Uw = tUw.ref();
261
262 if (primaryMesh().moving())
263 {
265
266 // Set up mapping values per patch
267
268 PtrMap<vectorField> patchValues(2*patches.size());
269
270 for (const label patchi : patches)
271 {
273 (
274 primaryMesh().boundaryMesh()[patchi]
275 );
276
277 if (wpp)
278 {
279 patchValues.set(patchi, wpp->Uwall());
280 }
281 }
282
283 if (patchValues.size())
284 {
285 // Map Up to area
286 tmp<vectorField> UsWall = vsmPtr_->mapToSurface(patchValues);
287
288 const vectorField& nHat =
290
291 Uw.primitiveFieldRef() = UsWall() - nHat*(UsWall() & nHat);
293 }
294
295 return tUw;
296}
297
298
300{
301 auto tUs = areaVectorField::New
302 (
303 "tUs",
305 regionMesh(),
307 );
308
309 // Uf quadratic profile
310 tUs.ref() = Foam::sqrt(2.0)*Uf_;
311
312 return tUs;
313}
314
315
317{
318 const volVectorField& Uprimary =
320
321 auto tUp = areaVectorField::New
322 (
323 "tUp",
325 regionMesh(),
327 );
328 auto& Up = tUp.ref();
329
330
331 // Set up mapping values per patch
332
334
335 PtrMap<vectorField> patchValues(2*patches.size());
336
337 // U tangential on primary
338 for (const label patchi : patches)
339 {
340 const fvPatchVectorField& Uw = Uprimary.boundaryField()[patchi];
341
342 patchValues.set(patchi, -Uw.snGrad());
343 }
344
345
346 // Map U tang to surface
347 vsmPtr_->mapToSurface(patchValues, Up.primitiveFieldRef());
348
349 Up.primitiveFieldRef() *= h_.primitiveField();
350
351 // U tangent on surface
354 Up.primitiveFieldRef() -= nHat*(Up.primitiveField() & nHat);
355
356 return tUp;
357}
358
359
361{
362 auto tpg = areaScalarField::New
363 (
364 "tpg",
366 regionMesh(),
368 );
369 auto& pfg = tpg.ref();
370
371 if (!pName_.empty())
372 {
373 vsmPtr_->mapInternalToSurface
374 (
375 primaryMesh().lookupObject<volScalarField>(pName_),
376 pfg.primitiveFieldRef()
377 );
378 }
379
380 return tpg;
381}
382
383
385{
386 auto talpha = areaScalarField::New
387 (
388 "talpha",
390 regionMesh(),
392 );
393 auto& alpha = talpha.ref();
395 alpha = pos0(h_ - deltaWet_);
396
397 return talpha;
398}
399
400
402(
403 const label patchi,
404 const label facei,
405 const scalar massSource,
406 const vector& momentumSource,
407 const scalar pressureSource,
408 const scalar energySource
409)
411 massSource_.boundaryFieldRef()[patchi][facei] += massSource;
412 pnSource_.boundaryFieldRef()[patchi][facei] += pressureSource;
413 momentumSource_.boundaryFieldRef()[patchi][facei] += momentumSource;
414}
415
418{
420}
421
422
424{
425 if (debug && primaryMesh().time().writeTime())
426 {
430 }
431
435
437}
438
439
440// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
441
442} // End namespace areaSurfaceFilmModels
443} // End namespace regionModels
444} // End namespace Foam
445
446// ************************************************************************* //
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())
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
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.
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
A HashTable of pointers to objects of type <T> with a label key.
Definition PtrMap.H:49
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
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 dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition dictionary.H:487
const edgeScalarField & deltaCoeffs() const
Return reference to difference factors array.
const labelList & whichPolyPatches() const
The polyPatches related to the areaMesh, in sorted order.
Definition faMeshI.H:128
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition faMesh.C:1189
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
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 word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition polyMesh.C:832
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static autoPtr< liquidFilmBase > New(const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected model using dictionary.
virtual void addSources(const label patchi, const label facei, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource=0)
Add sources.
tmp< areaScalarField > alpha() const
Wet indicator using h0.
label nFilmCorr_
Number of film thickness correctors.
liquidFilmBase(const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from type name and mesh and dict.
Switch momentumPredictor_
Flag to enable momentum predictor.
tmp< areaScalarField > pg() const
Map primary static pressure.
volScalarField pnSource_
Normal pressure by particles.
dimensionedScalar h0_
Smallest film thickness.
tmp< areaVectorField > Uw() const
Wall velocity.
tmp< areaVectorField > Up() const
Primary region velocity at film hight. Assume the film to be.
label nCorr_
Number of PISO-like inner correctors.
dimensionedScalar deltaWet_
Film thickness beyond which face is assumed to be wet.
virtual scalar CourantNumber() const
Courant number evaluation.
tmp< areaVectorField > Us() const
Film surface film velocity.
virtual void postEvolveRegion()
Post-evolve region.
const dictionary & solution() const
Return the solution dictionary.
autoPtr< volSurfaceMapping > vsmPtr_
Volume/surface mapping.
const word & areaName() const noexcept
The finite-area mesh name (extracted from dictionary).
const Time & time() const noexcept
Return the reference to the time database.
virtual void preEvolveRegion()
Pre-evolve region.
const faMesh & regionMesh() const
Return the region mesh database.
const fvMesh & primaryMesh() const noexcept
Return the reference to the primary mesh database.
regionFaModel(const fvMesh &mesh, const word &regionType, const word &modelName, const dictionary &dict, bool readFields=true)
Construct from mesh and name and dict.
A class for managing temporary objects.
Definition tmp.H:75
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
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
word timeName
Definition getTimeIndex.H:3
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for finite-area.
Definition limitHeight.C:30
const Foam::word liquidFilmName("liquidFilm")
Namespace for OpenFOAM.
const dimensionSet dimPressure
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< vector, faPatchField, areaMesh > areaVectorField
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimAcceleration
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
fvPatchField< vector > fvPatchVectorField
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Calculate the second temporal derivative.
volScalarField & alpha
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & e
scalar CoNum
scalar velMag