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-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 "liquidFilmBase.H"
29#include "faMesh.H"
30#include "gravityMeshObject.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace regionModels
41{
42namespace areaSurfaceFilmModels
43{
44
45// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46
48
50
51const Foam::word liquidFilmName("liquidFilm");
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
56(
57 const word& modelType,
58 const fvPatch& p,
59 const dictionary& dict
60)
61:
62 regionFaModel(p, liquidFilmName, modelType, dict, true),
63
64 momentumPredictor_
65 (
66 this->solution().subDict("PIMPLE").get<bool>("momentumPredictor")
67 ),
68 nOuterCorr_
69 (
70 this->solution().subDict("PIMPLE").get<label>("nOuterCorr")
71 ),
72 nCorr_(this->solution().subDict("PIMPLE").get<label>("nCorr")),
73 nFilmCorr_
74 (
75 this->solution().subDict("PIMPLE").get<label>("nFilmCorr")
76 ),
77
78 h0_("h0", dimLength, 1e-7, dict),
79
80 deltaWet_("deltaWet", dimLength, 1e-4, dict),
81
82 UName_(dict.get<word>("U")),
83
84 pName_(dict.lookupOrDefault<word>("p", word::null)),
85
86 pRef_(dict.get<scalar>("pRef")),
87
88 h_
89 (
91 (
92 "hf_" + regionName_,
93 primaryMesh().time().timeName(),
94 primaryMesh(),
95 IOobject::MUST_READ,
96 IOobject::AUTO_WRITE
97 ),
98 regionMesh()
99 ),
100
101 Uf_
102 (
104 (
105 "Uf_" + regionName_,
106 primaryMesh().time().timeName(),
107 primaryMesh(),
108 IOobject::MUST_READ,
109 IOobject::AUTO_WRITE
110 ),
111 regionMesh()
112 ),
113 pf_
114 (
116 (
117 "pf_" + regionName_,
118 primaryMesh().time().timeName(),
119 primaryMesh(),
120 IOobject::NO_READ,
121 IOobject::AUTO_WRITE
122 ),
123 regionMesh(),
125 ),
126 ppf_
127 (
129 (
130 "ppf_" + regionName_,
131 primaryMesh().time().timeName(),
132 primaryMesh(),
133 IOobject::NO_READ,
134 IOobject::NO_WRITE
135 ),
136 regionMesh(),
138 ),
139 phif_
140 (
142 (
143 "phif_" + regionName_,
144 primaryMesh().time().timeName(),
145 primaryMesh(),
146 IOobject::READ_IF_PRESENT,
147 IOobject::AUTO_WRITE
148 ),
149 fac::interpolate(Uf_) & regionMesh().Le()
150 ),
151
152 phi2s_
153 (
155 (
156 "phi2s_" + regionName_,
157 primaryMesh().time().timeName(),
158 primaryMesh(),
159 IOobject::READ_IF_PRESENT,
160 IOobject::AUTO_WRITE
161 ),
162 fac::interpolate(h_*Uf_) & regionMesh().Le()
163 ),
164
165
166 gn_
167 (
169 (
170 "gn",
171 primaryMesh().time().timeName(),
172 primaryMesh(),
173 IOobject::NO_READ,
174 IOobject::NO_WRITE
175 ),
176 regionMesh(),
178 ),
179
180 g_(meshObjects::gravity::New(primaryMesh().time())),
181
182 massSource_
183 (
185 (
186 "massSource",
187 primaryMesh().time().timeName(),
188 primaryMesh()
189 ),
190 primaryMesh(),
192 calculatedFvPatchField<scalar>::typeName
193 ),
194
195 momentumSource_
196 (
198 (
199 "momentumSource",
200 primaryMesh().time().timeName(),
201 primaryMesh()
202 ),
203 primaryMesh(),
206 ),
207
208 pnSource_
209 (
211 (
212 "pnSource",
213 primaryMesh().time().timeName(),
214 primaryMesh()
215 ),
216 primaryMesh(),
218 calculatedFvPatchField<scalar>::typeName
219 ),
220
221 energySource_
222 (
224 (
225 "energySource",
226 primaryMesh().time().timeName(),
227 primaryMesh()
228 ),
229 primaryMesh(),
231 calculatedFvPatchField<scalar>::typeName
232 ),
233
234 addedMassTotal_(0),
235
236 faOptions_(Foam::fa::options::New(p))
237{
239
240 gn_ = g_ & ns;
241
242 if (!faOptions_.optionList::size())
243 {
244 Info << "No finite area options present" << endl;
245 }
246}
247
248
249// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
250
252{}
253
254
255// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
256
258{
259 scalar CoNum = 0.0;
260 scalar velMag = 0.0;
261
262 edgeScalarField SfUfbyDelta
263 (
265 );
266
267 CoNum =
268 max(SfUfbyDelta/regionMesh().magLe()).value()*time().deltaT().value();
269
270 velMag = max(mag(phif_)/regionMesh().magLe()).value();
271
274
275 Info<< "Max film Courant Number: " << CoNum
276 << " Film velocity magnitude: " << velMag << endl;
277
278 return CoNum;
279}
280
281
283{
285 (
287 (
289 (
290 "tUw",
293 ),
294 regionMesh(),
296 )
297 );
298
299 areaVectorField& Uw = tUw.ref();
300
301 const polyPatch& pp = primaryMesh().boundaryMesh()[patch_.index()];
302
303 if
304 (
305 primaryMesh().moving()
306 && isA<movingWallVelocityFvPatchVectorField>(pp)
307 )
308 {
310 refCast<const movingWallVelocityFvPatchVectorField>(pp);
311
312 tmp<vectorField> tUwall = wpp.Uwall();
313
314 // Map Up to area
315 tmp<vectorField> UsWall = vsmPtr_->mapToSurface(tUwall());
316
317 const vectorField& nHat =
319
320 Uw.primitiveFieldRef() = UsWall() - nHat*(UsWall() & nHat);
321 }
322
323 return tUw;
324}
325
326
328{
330 (
332 (
334 (
335 "tUs",
338 ),
339 regionMesh(),
341 )
342 );
343
344 // Uf quadratic profile
345 tUs.ref() = Foam::sqrt(2.0)*Uf_;
346
347 return tUs;
348}
349
350
352{
353 const label patchi = patch_.index();
354
355 const volVectorField& Uprimary =
357
358 const fvPatchVectorField& Uw = Uprimary.boundaryField()[patchi];
359
361 (
363 (
365 (
366 "tUp",
369 ),
370 regionMesh(),
372 )
373 );
374
375 areaVectorField& Up = tUp.ref();
376
378
379 // map areas h to hp on primary
380 vsmPtr_->mapToField(h_, hp);
381
383
384 // U tangential on primary
385 const vectorField Ust(-Uw.snGrad()*hp);
386
387 // Map U tang to surface
388 Up.primitiveFieldRef() = vsmPtr_->mapToSurface(Ust);
389
390 // U tangent on surface
391 Up.primitiveFieldRef() -= nHat*(Up.primitiveField() & nHat);
392
393 return tUp;
394}
395
396
398{
400 (
402 (
404 (
405 "tpg",
408 ),
409 regionMesh(),
411 )
412 );
413
414 areaScalarField& pfg = tpg.ref();
415
416 if (pName_ != word::null)
417 {
418 const volScalarField& pp =
420
422 const_cast<volScalarField::Boundary&>(pp.boundaryField());
423
424 //pw -= pRef_;
425
426 pfg.primitiveFieldRef() = vsmPtr_->mapInternalToSurface<scalar>(pw)();
427 }
428 return tpg;
429}
430
431
433{
435 (
437 (
439 (
440 "talpha",
443 ),
444 regionMesh(),
446 )
447 );
448 areaScalarField& alpha = talpha.ref();
449
450 alpha = pos0(h_ - deltaWet_);
451
452 return talpha;
453}
454
455
457(
458 const label patchi,
459 const label facei,
460 const scalar massSource,
461 const vector& momentumSource,
462 const scalar pressureSource,
463 const scalar energySource
464)
465{
466 massSource_.boundaryFieldRef()[patchi][facei] += massSource;
467 pnSource_.boundaryFieldRef()[patchi][facei] += pressureSource;
468 momentumSource_.boundaryFieldRef()[patchi][facei] += momentumSource;
469}
470
471
473{
475}
476
477
479{
480 if (debug && primaryMesh().time().writeTime())
481 {
485 }
486
490
492}
493
494
496{
497 return faOptions_;
498}
499
500
502{
503 return Uf_;
504}
505
506
508{
509 return gn_;
510}
511
512
514{
515 return g_;
516}
517
518
520{
521 return h_;
522}
523
524
526{
527 return phif_;
528}
529
530
532{
533 return phi2s_;
534}
535
536
538{
539 return h0_;
540}
541
542
544{
545 return *this;
546}
547
548
550{
551 return pRef_;
552}
553
554
556{
557 return UName_;
558}
559
560
561// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
562
563} // End namespace areaSurfaceFilmModels
564} // End namespace regionModels
565} // End namespace Foam
566
567// ************************************************************************* //
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary 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
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:55
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const Type & value() const
Return const reference to value.
const edgeScalarField & deltaCoeffs() const
Return reference to difference factors array.
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition: faMesh.C:830
Finite-area options.
Definition: faOptions.H:58
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual label size() const
Return size.
Definition: fvPatch.H:185
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:203
This boundary condition provides a velocity condition for cases with moving walls.
tmp< vectorField > Uwall() const
Return wall velocity field.
const Type & lookupObject(const word &name, const bool recursive=false) const
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
virtual bool write(const bool valid=true) const
Write using setting from DB.
const edgeScalarField & phif() const
Access to momentum flux.
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.
const regionFaModel & region() const
Access to this region.
Foam::fa::options & faOptions()
Return faOptions.
tmp< areaScalarField > pg() const
Map primary static pressure.
volScalarField pnSource_
Normal pressure by particles.
dimensionedScalar h0_
Smallest numerical thickness.
tmp< areaVectorField > Uw() const
Wall velocity.
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.
const edgeScalarField & phi2s() const
Access continuity flux.
const dimensionedScalar & h0() const
Return h0.
dimensionedScalar deltaWet_
Delta wet for sub-models.
virtual scalar CourantNumber() const
Courant number evaluation.
tmp< areaVectorField > Us() const
Film surface film velocity.
const uniformDimensionedVectorField & g() const
Gravity.
const areaScalarField & gn() const
Access const reference gn.
Base class for area region models.
virtual void postEvolveRegion()
Post-evolve region.
const Time & time() const
Return the reference to the time database.
autoPtr< volSurfaceMapping > vsmPtr_
Volume-to surface mapping.
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
virtual void preEvolveRegion()
Pre-evolve region.
const fvPatch & patch_
Reference to fvPatch.
const faMesh & regionMesh() const
Return the region mesh database.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:66
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
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
volScalarField & p
bool
Definition: EEqn.H:20
word timeName
Definition: getTimeIndex.H:3
const Foam::word liquidFilmName("liquidFilm")
Namespace for OpenFOAM.
const dimensionSet dimPressure
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
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:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const dimensionSet dimAcceleration
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Calculate the second temporal derivative.
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & e
Definition: createFields.H:11
scalar CoNum
scalar velMag