Loading...
Searching...
No Matches
OwenRyleyModel.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2021-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "OwenRyleyModel.H"
30#include "processorFaPatch.H"
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
38{
41
42
43// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44
45tmp<areaScalarField> OwenRyleyModel::calcInvR1
46(
47 const areaVectorField& U
48) const
49{
50 const dimensionedScalar smallU(dimVelocity, ROOTVSMALL);
51 const areaVectorField UHat(U/(mag(U) + smallU));
52
53 auto tinvR1 = areaScalarField::New
54 (
55 "invR1",
57 (UHat & (UHat & -gradNHat_))
58 );
59 scalarField& invR1 = tinvR1.ref().primitiveFieldRef();
60
61 // Apply defined patch radii
62 if (definedPatchRadii_ > 0)
63 {
64 invR1 = 1.0/max(1e-6, definedPatchRadii_);
65 }
66
67 // Filter out large radii
68 for (auto& x : invR1)
69 {
70 if (mag(x) < 1e-6)
71 {
72 x = -1;
73 }
74 }
75
76 return tinvR1;
77}
78
79
80tmp<scalarField> OwenRyleyModel::calcCosAngle
81(
82 const edgeScalarField& phi,
83 const scalarField& invR1
84) const
85{
86 const areaVectorField& U = film().Uf();
87 const dimensionedScalar smallU(dimVelocity, ROOTVSMALL);
88 const areaVectorField UHat(U/(mag(U) + smallU));
89
90 const vector gHat(normalised(film().g().value()));
91
92 const faMesh& mesh = film().regionMesh();
93 const labelUList& own = mesh.edgeOwner();
94 const labelUList& nbr = mesh.edgeNeighbour();
95
96 scalarField phiMax(mesh.nFaces(), -GREAT);
97 scalarField cosAngle(UHat.size(), Zero);
98
99 // Internal edges
100 forAll(nbr, edgei)
101 {
102 const label cellO = own[edgei];
103 const label cellN = nbr[edgei];
104 const scalar phiEdge = phi[edgei];
105
106 if (phiEdge > phiMax[cellO])
107 {
108 phiMax[cellO] = phiEdge;
109 cosAngle[cellO] = -gHat & UHat[cellN];
110 }
111 if (-phiEdge > phiMax[cellN])
112 {
113 phiMax[cellN] = -phiEdge;
114 cosAngle[cellN] = -gHat & UHat[cellO];
115 }
116 }
117
118 // Processor edges
119 for (const auto& phip : phi.boundaryField())
120 {
121 if (isA<processorFaPatch>(phip.patch()))
122 {
123 const auto& edgeFaces = phip.patch().edgeFaces();
124 const auto& UHatp = UHat.boundaryField()[phip.patch().index()];
125 forAll(phip, edgei)
126 {
127 const label cellO = edgeFaces[edgei];
128 if (phip[edgei] > phiMax[cellO])
129 {
130 phiMax[cellO] = phip[edgei];
131 cosAngle[cellO] = -gHat & UHatp[edgei];
132 }
133 }
134 }
135 }
136
137 cosAngle *= pos(invR1);
138
139
140 if (debug && mesh.time().writeTime())
141 {
142 {
143 areaScalarField areaCosAngle
144 (
145 mesh.newIOobject("cosAngle"),
146 mesh,
148 );
149 areaCosAngle.primitiveFieldRef() = cosAngle;
150 areaCosAngle.correctBoundaryConditions();
151 areaCosAngle.write();
152 }
153
154 {
155 areaScalarField areaInvR1
156 (
157 mesh.newIOobject("InvR1"),
158 mesh,
160 );
161 areaInvR1.primitiveFieldRef() = invR1;
162 areaInvR1.write();
163 }
164 }
165
166 return clamp(cosAngle, scalarMinMax(-1, 1));
167}
168
169
170tmp<scalarField> OwenRyleyModel::netForce() const
171{
172 const faMesh& mesh = film().regionMesh();
173
174 const areaScalarField& h = film().h();
175 const areaVectorField& U = film().Uf();
176 const edgeScalarField& phi = film().phi2s();
177 const areaScalarField& rho = film().rho();
178 const areaScalarField& sigma = film().sigma();
179
180 const scalarField magSqrU(magSqr(U));
181 const scalarField invR1(calcInvR1(U));
182 const scalarField cosAngle(calcCosAngle(phi, invR1));
183
184 const scalar magG = mag(film().g().value());
185
186 // Initialize the net-force magnitude
187 auto tFnet = tmp<scalarField>::New(mesh.nFaces(), Zero);
188 auto& Fnet = tFnet.ref();
189
190 forAll(Fnet, i)
191 {
192 if ((invR1[i] > minInvR1()) && (h[i]*invR1[i] > minHbyR1()))
193 {
194 const scalar R1 = 1.0/(invR1[i] + ROOTVSMALL);
195 const scalar R2 = R1 + h[i];
196
197 // Inertial force (OR:Eq. 11; '2' is an exponent in the Eq.)
198 const scalar Fi = -4.0/3.0*h[i]*rho[i]*magSqrU[i]*invR1[i];
199
200 // Body force (OR:Eq. 11)
201 const scalar Fb =
202 -0.5*rho[i]*magG*invR1[i]*(sqr(R1) - sqr(R2))*cosAngle[i];
203
204 // Surface force (OR:Eq. 11)
205 const scalar Fs = sigma[i]/R2;
206
207 Fnet[i] = Fi + Fb + Fs;
208 }
209 }
210
211 if (debug && mesh.time().writeTime())
212 {
213 {
214 areaScalarField areaFnet
215 (
216 mesh.newIOobject("Fnet"),
217 mesh,
219 );
220 areaFnet.primitiveFieldRef() = Fnet;
221 areaFnet.write();
222 }
223 }
225 return tFnet;
226}
227
228
229// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
230
232(
234 const dictionary& dict
235)
236:
238 fThreshold_(dict.getOrDefault<scalar>("fThreshold", 1e-8)),
239 definedPatchRadii_(dict.getOrDefault<scalar>("definedPatchRadii", 0)),
240 minHbyR1_(dict.getOrDefault<scalar>("deltaByR1Min", 0)),
241 minInvR1_(dict.getOrDefault<scalar>("minInvR1", 5)),
242 gradNHat_(fac::grad(film.regionMesh().faceAreaNormals()))
243{
244 if (mag(film.g().value()) < ROOTVSMALL)
245 {
246 FatalErrorInFunction
247 << "Gravitational acceleration must be non-zero."
248 << exit(FatalError);
249 }
250}
251
252
253// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
254
256{
257 const faMesh& mesh = film().regionMesh();
258
259 tmp<scalarField> tFnet = netForce();
260 const auto& Fnet = tFnet.cref();
261
262 // Initialize the mass ratio of film separation
263 auto tseparated = tmp<scalarField>::New(mesh.nFaces(), Zero);
264 auto& separated = tseparated.ref();
265
266 forAll(Fnet, i)
267 {
268 if ((Fnet[i] + fThreshold_) < 0)
269 {
270 separated[i] = 1;
271 }
272 }
273
274 return tseparated;
275}
276
277
278// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279
280} // End namespace filmSeparationModels
281} // End namespace Foam
282
283
284// ************************************************************************* //
285
CGAL::indexedFace< K > Fb
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
static tmp< GeometricField< scalar, faPatchField, areaMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=faPatchField< scalar >::calculatedType())
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
void correctBoundaryConditions()
Correct boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const Type & value() const noexcept
Return const reference to value.
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition faMesh.H:140
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 round edges for full separation (Owen & Ryley,...
OwenRyleyModel(const regionModels::areaSurfaceFilmModels::liquidFilmBase &film, const dictionary &dict)
Construct from the base film model and dictionary.
virtual tmp< scalarField > separatedMassRatio() const
Calculate the mass ratio of film separation.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
virtual const areaScalarField & sigma() const =0
Access const reference sigma.
const uniformDimensionedVectorField & g() const noexcept
Gravity.
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.
const faMesh & regionMesh() const
Return the region mesh database.
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
U
Definition pEqn.H:72
dynamicFvMesh & mesh
Namespace for handling debugging switches.
Definition debug.C:45
A namespace for various filmSeparation model implementations.
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimVelocity
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
GeometricField< vector, faPatchField, areaMesh > areaVectorField
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
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
const dimensionSet dimForce
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
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
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Calculate the second temporal derivative.
dictionary dict
volScalarField & h
volScalarField & e
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299