Loading...
Searching...
No Matches
FriedrichModel.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) 2024-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 "FriedrichModel.H"
29#include "processorFaPatch.H"
30#include "unitConversion.H"
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
38{
41
42
43const Foam::Enum
44<
45 FriedrichModel::separationType
46>
47FriedrichModel::separationTypeNames
48({
49 { separationType::FULL, "full" },
50 { separationType::PARTIAL , "partial" },
51});
52
53
54// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55
56bitSet FriedrichModel::calcCornerEdges() const
57{
58 bitSet cornerEdges(mesh().nEdges(), false);
59
60 const areaVectorField& faceCentres = mesh().areaCentres();
61 const areaVectorField& faceNormals = mesh().faceAreaNormals();
62
63 const labelUList& own = mesh().edgeOwner();
64 const labelUList& nbr = mesh().edgeNeighbour();
65
66 // Check if internal face-normal vectors diverge (no separation)
67 // or converge (separation may occur)
68 forAll(nbr, edgei)
69 {
70 const label faceO = own[edgei];
71 const label faceN = nbr[edgei];
72
73 cornerEdges[edgei] = isCornerEdgeSharp
74 (
75 faceCentres[faceO],
76 faceCentres[faceN],
77 faceNormals[faceO],
78 faceNormals[faceN]
79 );
80 }
81
82
83 // Skip the rest of the routine if the simulation is a serial run
84 if (!Pstream::parRun()) return cornerEdges;
85
86 // Check if processor face-normal vectors diverge (no separation)
87 // or converge (separation may occur)
89
90 for (const faPatch& fap : patches)
91 {
92 if (isA<processorFaPatch>(fap))
93 {
94 const label patchi = fap.index();
95 const auto& edgeFaces = fap.edgeFaces();
96 const label internalEdgei = fap.start();
97
98 const auto& faceCentresp = faceCentres.boundaryField()[patchi];
99 const auto& faceNormalsp = faceNormals.boundaryField()[patchi];
100
101 forAll(faceNormalsp, bndEdgei)
102 {
103 const label faceO = edgeFaces[bndEdgei];
104 const label meshEdgei = internalEdgei + bndEdgei;
105
106 cornerEdges[meshEdgei] = isCornerEdgeSharp
107 (
108 faceCentres[faceO],
109 faceCentresp[bndEdgei],
110 faceNormals[faceO],
111 faceNormalsp[bndEdgei]
112 );
113 }
114 }
115 }
116
117 return cornerEdges;
118}
119
120
121bool FriedrichModel::isCornerEdgeSharp
122(
123 const vector& faceCentreO,
124 const vector& faceCentreN,
125 const vector& faceNormalO,
126 const vector& faceNormalN
127) const
128{
129 // Calculate the relative position of centres of faces sharing an edge
130 const vector relativePosition(faceCentreN - faceCentreO);
131
132 // Calculate the relative normal of faces sharing an edge
133 const vector relativeNormal(faceNormalN - faceNormalO);
134
135 // Return true if the face normals converge, meaning that the edge is sharp
136 return ((relativeNormal & relativePosition) < -1e-8);
137}
138
139
140scalarList FriedrichModel::calcCornerAngles() const
141{
142 scalarList cornerAngles(mesh().nEdges(), Zero);
143
144 const areaVectorField& faceNormals = mesh().faceAreaNormals();
145
146 const labelUList& own = mesh().edgeOwner();
147 const labelUList& nbr = mesh().edgeNeighbour();
148
149 // Process internal edges
150 forAll(nbr, edgei)
151 {
152 if (!cornerEdges_[edgei]) continue;
153
154 const label faceO = own[edgei];
155 const label faceN = nbr[edgei];
156
157 cornerAngles[edgei] = calcCornerAngle
158 (
159 faceNormals[faceO],
160 faceNormals[faceN]
161 );
162 }
163
164
165 // Skip the rest of the routine if the simulation is a serial run
166 if (!Pstream::parRun()) return cornerAngles;
167
168 // Process processor edges
169 const faBoundaryMesh& patches = mesh().boundary();
170
171 for (const faPatch& fap : patches)
172 {
173 if (isA<processorFaPatch>(fap))
174 {
175 const label patchi = fap.index();
176 const auto& edgeFaces = fap.edgeFaces();
177 const label internalEdgei = fap.start();
178
179 const auto& faceNormalsp = faceNormals.boundaryField()[patchi];
180
181 forAll(faceNormalsp, bndEdgei)
182 {
183 const label faceO = edgeFaces[bndEdgei];
184 const label meshEdgei = internalEdgei + bndEdgei;
185
186 if (!cornerEdges_[meshEdgei]) continue;
187
188 cornerAngles[meshEdgei] = calcCornerAngle
189 (
190 faceNormals[faceO],
191 faceNormalsp[bndEdgei]
192 );
193 }
194 }
195 }
196
197 return cornerAngles;
198}
199
200
201scalar FriedrichModel::calcCornerAngle
202(
203 const vector& faceNormalO,
204 const vector& faceNormalN
205) const
206{
207 const scalar magFaceNormal = mag(faceNormalO)*mag(faceNormalN);
208
209 // Avoid any potential exceptions during the cosine calculations
210 if (magFaceNormal < SMALL) return 0;
211
212 scalar cosAngle = (faceNormalO & faceNormalN)/magFaceNormal;
213 cosAngle = clamp(cosAngle, -1, 1);
214
215 return std::acos(cosAngle);
216}
217
218
219bitSet FriedrichModel::calcSeparationFaces() const
220{
221 bitSet separationFaces(mesh().faces().size(), false);
222
223 const edgeScalarField& phis = film().phi2s();
224
225 const labelUList& own = mesh().edgeOwner();
226 const labelUList& nbr = mesh().edgeNeighbour();
227
228 // Process internal faces
229 forAll(nbr, edgei)
230 {
231 if (!cornerEdges_[edgei]) continue;
232
233 const label faceO = own[edgei];
234 const label faceN = nbr[edgei];
235
236 isSeparationFace
237 (
238 separationFaces,
239 phis[edgei],
240 faceO,
241 faceN
242 );
243 }
244
245
246 // Skip the rest of the routine if the simulation is a serial run
247 if (!Pstream::parRun()) return separationFaces;
248
249 // Process processor faces
250 const faBoundaryMesh& patches = mesh().boundary();
251
252 for (const faPatch& fap : patches)
253 {
254 if (isA<processorFaPatch>(fap))
255 {
256 const label patchi = fap.index();
257 const auto& edgeFaces = fap.edgeFaces();
258 const label internalEdgei = fap.start();
259
260 const auto& phisp = phis.boundaryField()[patchi];
261
262 forAll(phisp, bndEdgei)
263 {
264 const label faceO = edgeFaces[bndEdgei];
265 const label meshEdgei(internalEdgei + bndEdgei);
266
267 if (!cornerEdges_[meshEdgei]) continue;
268
269 isSeparationFace
270 (
271 separationFaces,
272 phisp[bndEdgei],
273 faceO
274 );
275 }
276 }
277 }
278
279 return separationFaces;
280}
281
282
283void FriedrichModel::isSeparationFace
284(
285 bitSet& separationFaces,
286 const scalar phiEdge,
287 const label faceO,
288 const label faceN
289) const
290{
291 const scalar tol = 1e-8;
292
293 // Assuming there are no sources/sinks at the edge
294 if (phiEdge > tol) // From owner to neighbour
295 {
296 separationFaces[faceO] = true;
297 }
298 else if ((phiEdge < -tol) && (faceN != -1)) // From neighbour to owner
299 {
300 separationFaces[faceN] = true;
301 }
302}
303
304
305scalarList FriedrichModel::calcSeparationAngles
306(
307 const bitSet& separationFaces
308) const
309{
310 scalarList separationAngles(mesh().faces().size(), Zero);
311
312 const labelUList& own = mesh().edgeOwner();
313 const labelUList& nbr = mesh().edgeNeighbour();
314
315 // Process internal faces
316 forAll(nbr, edgei)
317 {
318 if (!cornerEdges_[edgei]) continue;
319
320 const label faceO = own[edgei];
321 const label faceN = nbr[edgei];
322
323 if (separationFaces[faceO])
324 {
325 separationAngles[faceO] = cornerAngles_[edgei];
326 }
327
328 if (separationFaces[faceN])
329 {
330 separationAngles[faceN] = cornerAngles_[edgei];
331 }
332 }
333
334
335 // Skip the rest of the routine if the simulation is a serial run
336 if (!Pstream::parRun()) return separationAngles;
337
338 // Process processor faces
339 const edgeScalarField& phis = film().phi2s();
340 const faBoundaryMesh& patches = mesh().boundary();
341
342 for (const faPatch& fap : patches)
343 {
344 if (isA<processorFaPatch>(fap))
345 {
346 const label patchi = fap.index();
347 const auto& edgeFaces = fap.edgeFaces();
348 const label internalEdgei = fap.start();
349
350 const auto& phisp = phis.boundaryField()[patchi];
351
352 forAll(phisp, bndEdgei)
353 {
354 const label faceO = edgeFaces[bndEdgei];
355 const label meshEdgei(internalEdgei + bndEdgei);
356
357 if (!cornerEdges_[meshEdgei]) continue;
358
359 if (separationFaces[faceO])
360 {
361 separationAngles[faceO] = cornerAngles_[meshEdgei];
362 }
363 }
364 }
365 }
366
367 return separationAngles;
368}
369
370
371tmp<scalarField> FriedrichModel::Fratio() const
372{
373 const areaVectorField Up(film().Up());
374 const areaVectorField& Uf = film().Uf();
375 const areaScalarField& h = film().h();
376 const areaScalarField& rho = film().rho();
377 const areaScalarField& mu = film().mu();
378 const areaScalarField& sigma = film().sigma();
379
380 // Identify the faces where separation may occur
381 const bitSet separationFaces(calcSeparationFaces());
382
383 // Calculate the corner angles corresponding to the separation faces
384 const scalarList separationAngles(calcSeparationAngles(separationFaces));
385
386 // Initialize the force ratio
387 auto tFratio = tmp<scalarField>::New(mesh().faces().size(), Zero);
388 auto& Fratio = tFratio.ref();
389
390 // Process internal faces
391 forAll(separationFaces, i)
392 {
393 // Skip the routine if the face is not a candidate for separation
394 if (!separationFaces[i]) continue;
395
396 // Calculate the corner-angle trigonometric values
397 const scalar sinAngle = std::sin(separationAngles[i]);
398 const scalar cosAngle = std::cos(separationAngles[i]);
399
400 // Reynolds number (FLW:Eq. 16)
401 const scalar Re = h[i]*mag(Uf[i])*rho[i]/mu[i];
402
403 // Weber number (FLW:Eq. 17)
404 const vector Urel(Up[i] - Uf[i]);
405 const scalar We = h[i]*rhop_*sqr(mag(Urel))/(2.0*sigma[i]);
406
407 // Characteristic breakup length (FLW:Eq. 15)
408 const scalar Lb =
409 0.0388*Foam::sqrt(h[i])*Foam::pow(Re, 0.6)*Foam::pow(We, -0.5);
410
411 // Force ratio - denominator (FLW:Eq. 20)
412 const scalar den =
413 sigma[i]*(sinAngle + 1.0) + rho[i]*magG_*h[i]*Lb*cosAngle;
414
415 if (mag(den) > 0)
416 {
417 // Force ratio (FLW:Eq. 20)
418 Fratio[i] = rho[i]*sqr(mag(Uf[i]))*h[i]*sinAngle/den;
419 }
420 }
421
422
423 // Skip the rest of the routine if the simulation is a serial run
424 if (!Pstream::parRun()) return tFratio;
425
426 // Process processor faces
427 const faBoundaryMesh& patches = mesh().boundary();
428
429 for (const faPatch& fap : patches)
430 {
431 if (isA<processorFaPatch>(fap))
432 {
433 const label patchi = fap.index();
434 const label internalEdgei = fap.start();
435
436 const auto& hp = h.boundaryField()[patchi];
437 const auto& Ufp = Uf.boundaryField()[patchi];
438 const auto& Upp = Up.boundaryField()[patchi];
439 const auto& rhop = rho.boundaryField()[patchi];
440 const auto& sigmap = sigma.boundaryField()[patchi];
441 const auto& mup = mu.boundaryField()[patchi];
442
443 forAll(hp, i)
444 {
445 // Skip the routine if the face is not a candidate for separation
446 if (!separationFaces[i]) continue;
447
448 const label meshEdgei = internalEdgei + i;
449
450 // Calculate the corner-angle trigonometric values
451 const scalar sinAngle = std::sin(cornerAngles_[meshEdgei]);
452 const scalar cosAngle = std::cos(cornerAngles_[meshEdgei]);
453
454 // Reynolds number (FLW:Eq. 16)
455 const scalar Re = hp[i]*mag(Ufp[i])*rhop[i]/mup[i];
456
457 // Weber number (FLW:Eq. 17)
458 const vector Urelp(Upp[i] - Ufp[i]);
459 const scalar We = hp[i]*rhop_*sqr(mag(Urelp))/(2.0*sigmap[i]);
460
461 // Characteristic breakup length (FLW:Eq. 15)
462 const scalar Lb =
463 0.0388*Foam::sqrt(hp[i])
464 *Foam::pow(Re, 0.6)*Foam::pow(We, -0.5);
465
466 // Force ratio - denominator (FLW:Eq. 20)
467 const scalar den =
468 sigmap[i]*(sinAngle + 1.0)
469 + rhop[i]*magG_*hp[i]*Lb*cosAngle;
470
471 if (mag(den) > 0)
472 {
473 // Force ratio (FLW:Eq. 20)
474 Fratio[i] = rhop[i]*sqr(mag(Ufp[i]))*hp[i]*sinAngle/den;
475 }
476 }
477 }
478 }
480 return tFratio;
481}
482
483
484// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
485
487(
489 const dictionary& dict
490)
491:
493 separation_
494 (
495 separationTypeNames.getOrDefault
496 (
497 "separationType",
498 dict,
499 separationType::FULL
500 )
501 ),
502 rhop_(dict.getScalar("rhop")),
503 magG_(mag(film.g().value())),
504 C0_(dict.getOrDefault<scalar>("C0", 0.882)),
505 C1_(dict.getOrDefault<scalar>("C1", -1.908)),
506 C2_(dict.getOrDefault<scalar>("C2", 1.264)),
507 cornerEdges_(calcCornerEdges()),
508 cornerAngles_(calcCornerAngles())
509{
510 if (rhop_ < VSMALL)
511 {
513 << "Primary-phase density, rhop: " << rhop_ << " must be non-zero."
515 }
516
517 if (mag(C2_) < VSMALL)
518 {
520 << "Empirical constant, C2 = " << C2_ << "cannot be zero."
522 }
523}
524
525
526// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
527
529{
530 tmp<scalarField> tFratio = Fratio();
531 const auto& Fratio = tFratio.cref();
532
533 // Initialize the mass ratio of film separation
534 auto tseparated = tmp<scalarField>::New(mesh().faces().size(), Zero);
535 auto& separated = tseparated.ref();
536
537
538 switch (separation_)
539 {
540 case separationType::FULL:
541 {
542 forAll(Fratio, i)
543 {
544 if (Fratio[i] > 1)
545 {
546 separated[i] = 1;
547 }
548 }
549 break;
550 }
551 case separationType::PARTIAL:
552 {
553 forAll(Fratio, i)
554 {
555 if (Fratio[i] > 1)
556 {
557 // (ZJD:Eq. 16)
558 separated[i] = C0_ + C1_*Foam::exp(-Fratio[i]/C2_);
559 }
560 }
561 break;
562 }
563 default:
564 break; // This should not happen.
565 }
566
567 if (debug && mesh().time().writeTime())
568 {
569 {
570 areaScalarField areaFratio
571 (
572 mesh().newIOobject("Fratio"),
573 mesh(),
575 );
576 areaFratio.primitiveFieldRef() = Fratio;
577 areaFratio.write();
578 }
579 }
580
581
582 return tseparated;
583}
584
585
586// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
587
588} // End namespace filmSeparationModels
589} // End namespace Foam
590
591
592// ************************************************************************* //
593
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const uniformDimensionedVectorField & g
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Finite area boundary mesh, which is a faPatch list with registered IO, a reference to the associated ...
const labelList & edgeNeighbour() const noexcept
Edge neighbour addressing.
Definition faMeshI.H:90
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition faMeshI.H:24
const labelList & edgeOwner() const noexcept
Edge owner addressing.
Definition faMeshI.H:84
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition faMesh.C:1189
const areaVectorField & areaCentres() const
Return face centres as areaVectorField.
Definition faMesh.C:1116
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition faPatch.H:76
A base class for filmSeparation models.
const regionModels::areaSurfaceFilmModels::liquidFilmBase & film() const
Return const access to the film properties.
const faMesh & mesh() const noexcept
Return const access to the finite-area mesh.
filmSeparationModel(const filmSeparationModel &)=delete
No copy construct.
Computes film-separation properties from sharp edges for full separation (Friedrich et al....
virtual tmp< scalarField > separatedMassRatio() const
Calculate the mass ratio of film separation.
FriedrichModel(const regionModels::areaSurfaceFilmModels::liquidFilmBase &film, const dictionary &dict)
Construct from the base film model and dictionary.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
virtual const areaScalarField & sigma() const =0
Access const reference sigma.
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.
const edgeScalarField & phi2s() const noexcept
Access continuity flux.
A class for managing temporary objects.
Definition tmp.H:75
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const polyBoundaryMesh & patches
autoPtr< surfaceVectorField > Uf
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
Urel
Definition pEqn.H:56
Namespace for handling debugging switches.
Definition debug.C:45
A namespace for various filmSeparation model implementations.
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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)
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
errorManip< error > abort(error &err)
Definition errorManip.H:139
const dimensionSet dimForce
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
scalarField Re(const UList< complex > &cmplx)
Extract real component.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
dictionary dict
volScalarField & h
volScalarField & e
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
edgeScalarField phis(IOobject("phis", runTime.timeName(), aMesh.thisDb(), IOobject::NO_READ, IOobject::NO_WRITE), linearEdgeInterpolate(Us) &aMesh.Le())
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.