curvatureSeparation.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) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "curvatureSeparation.H"
31 #include "Time.H"
32 #include "stringListOps.H"
33 #include "cyclicPolyPatch.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace regionModels
40 {
41 namespace areaSurfaceFilmModels
42 {
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 defineTypeNameAndDebug(curvatureSeparation, 0);
48 (
49  injectionModel,
52 );
53 
54 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
55 
57 (
58  const areaVectorField& U,
59  const scalarField& calcCosAngle
60 ) const
61 {
62  const dimensionedScalar smallU(dimVelocity, ROOTVSMALL);
63  const areaVectorField UHat(U/(mag(U) + smallU));
65  (
66  new areaScalarField("invR1", UHat & (UHat & -gradNHat_))
67  );
68 
69  scalarField& invR1 = tinvR1.ref().primitiveFieldRef();
70 
71  // apply defined patch radii
72  const scalar rMin = 1e-6;
73  const scalar definedInvR1 = 1.0/max(rMin, definedPatchRadii_);
74 
75  if (definedPatchRadii_ > 0)
76  {
77  invR1 = definedInvR1;
78  }
79 
80  // filter out large radii
81  const scalar rMax = 1e6;
82  forAll(invR1, i)
83  {
84  if ((mag(invR1[i]) < 1/rMax))
85  {
86  invR1[i] = -1.0;
87  }
88  }
89 
90  return tinvR1;
91 }
92 
93 
95 (
96  const edgeScalarField& phi
97 ) const
98 {
99  const areaVectorField& U = film().Uf();
100  const dimensionedScalar smallU(dimVelocity, ROOTVSMALL);
101  const areaVectorField UHat(U/(mag(U) + smallU));
102 
103  const faMesh& mesh = film().regionMesh();
104  const labelUList& own = mesh.edgeOwner();
105  const labelUList& nbr = mesh.edgeNeighbour();
106 
107  scalarField phiMax(mesh.nFaces(), -GREAT);
108  scalarField cosAngle(UHat.size(), Zero);
109 
110  const scalarField invR1(calcInvR1(U, cosAngle));
111 
112  forAll(nbr, edgei)
113  {
114  const label cellO = own[edgei];
115  const label cellN = nbr[edgei];
116 
117  if (phi[edgei] > phiMax[cellO])
118  {
119  phiMax[cellO] = phi[edgei];
120  cosAngle[cellO] = -gHat_ & UHat[cellN];
121  }
122  if (-phi[edgei] > phiMax[cellN])
123  {
124  phiMax[cellN] = -phi[edgei];
125  cosAngle[cellN] = -gHat_ & UHat[cellO];
126  }
127  }
128 
129  cosAngle *= pos(invR1);
130 
131  // checks
132  if (debug && mesh.time().writeTime())
133  {
134  areaScalarField volCosAngle
135  (
136  IOobject
137  (
138  "cosAngle",
139  mesh.time().timeName(),
140  mesh.thisDb(),
142  ),
143  mesh,
145  );
146  volCosAngle.primitiveFieldRef() = cosAngle;
147  volCosAngle.correctBoundaryConditions();
148  volCosAngle.write();
149  }
150 
151  return clamp(cosAngle, scalarMinMax(-1, 1));
152 }
153 
154 
155 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
156 
157 curvatureSeparation::curvatureSeparation
158 (
159  liquidFilmBase& film,
160  const dictionary& dict
161 )
162 :
163  injectionModel(type(), film, dict),
164  gradNHat_(fac::grad(film.regionMesh().faceAreaNormals())),
165  deltaByR1Min_(coeffDict_.getOrDefault<scalar>("deltaByR1Min", 0)),
166  definedPatchRadii_
167  (
168  coeffDict_.getOrDefault<scalar>("definedPatchRadii", 0)
169  ),
170  magG_(mag(film.g().value())),
171  gHat_(Zero),
172  fThreshold_
173  (
174  coeffDict_.getOrDefault<scalar>("fThreshold", 1e-8)
175  ),
176  minInvR1_
177  (
178  coeffDict_.getOrDefault<scalar>("minInvR1", 5)
179  )
180 {
181  if (magG_ < ROOTVSMALL)
182  {
184  << "Acceleration due to gravity must be non-zero"
185  << exit(FatalError);
186  }
187 
189 }
190 
191 
192 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
193 
195 (
196  scalarField& availableMass,
197  scalarField& massToInject,
198  scalarField& diameterToInject
199 )
200 {
201  const faMesh& mesh = film().regionMesh();
202 
203  const areaScalarField& delta = film().h();
204  const areaVectorField& U = film().Uf();
205  const edgeScalarField& phi = film().phi2s();
206  const areaScalarField& rho = film().rho();
207  const scalarField magSqrU(magSqr(film().Uf()));
208  const areaScalarField& sigma = film().sigma();
209 
210  const scalarField cosAngle(calcCosAngle(phi));
211  const scalarField invR1(calcInvR1(U, cosAngle));
212 
213 
214  // calculate force balance
215  scalarField Fnet(mesh.nFaces(), Zero);
216  scalarField separated(mesh.nFaces(), Zero);
217 
218  forAll(invR1, i)
219  {
220  if ((invR1[i] > minInvR1_) && (delta[i]*invR1[i] > deltaByR1Min_))
221  {
222  const scalar R1 = 1.0/(invR1[i] + ROOTVSMALL);
223  const scalar R2 = R1 + delta[i];
224 
225  // inertial force
226  const scalar Fi = -delta[i]*rho[i]*magSqrU[i]*72.0/60.0*invR1[i];
227 
228  // body force
229  const scalar Fb =
230  - 0.5*rho[i]*magG_*invR1[i]*(sqr(R1) - sqr(R2))*cosAngle[i];
231 
232  // surface force
233  const scalar Fs = sigma[i]/R2;
234 
235  Fnet[i] = Fi + Fb + Fs;
236 
237  if (Fnet[i] + fThreshold_ < 0)
238  {
239  separated[i] = 1.0;
240  }
241  }
242  }
243 
244  // inject all available mass
245  massToInject = separated*availableMass;
246  diameterToInject = separated*delta;
247  availableMass -= separated*availableMass;
248 
249  addToInjectedMass(sum(massToInject));
250 
251  if (debug && mesh.time().writeTime())
252  {
253  areaScalarField volFnet
254  (
255  IOobject
256  (
257  "Fnet",
258  mesh.time().timeName(),
259  mesh.thisDb(),
261  ),
262  mesh,
264  );
265  volFnet.primitiveFieldRef() = Fnet;
266  volFnet.write();
267 
268  areaScalarField areaSeparated
269  (
270  IOobject
271  (
272  "separated",
273  mesh.time().timeName(),
274  mesh.thisDb(),
276  ),
277  mesh,
279  );
280  areaSeparated.primitiveFieldRef() = separated;
281  areaSeparated.write();
282 
283  areaScalarField areaMassToInject
284  (
285  IOobject
286  (
287  "massToInject",
288  mesh.time().timeName(),
289  mesh.thisDb(),
291  ),
292  mesh,
294  );
295  areaMassToInject.primitiveFieldRef() = massToInject;
296  areaMassToInject.write();
297 
298  areaScalarField areaInvR1
299  (
300  IOobject
301  (
302  "InvR1",
303  mesh.time().timeName(),
304  mesh.thisDb(),
306  ),
307  mesh,
309  );
310  areaInvR1.primitiveFieldRef() = invR1;
311  areaInvR1.write();
312  }
313 
315 }
316 
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 } // End namespace areaSurfaceFilmModels
321 } // End namespace regionModels
322 } // End namespace Foam
323 
324 // ************************************************************************* //
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:87
scalar delta
const Type & value() const noexcept
Return const reference to value.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
void addToInjectedMass(const scalar dMass)
Add to injected mass.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
Calculate the second temporal derivative.
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition: MinMax.H:97
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const areaScalarField & h() const
Access const reference h.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
virtual const areaScalarField & sigma() const =0
Access const reference sigma.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
scalar deltaByR1Min_
Minimum gravity driven film thickness (non-dimensionalised delta/R1)
Operations on lists of strings.
const uniformDimensionedVectorField & g() const
Gravity.
bool writeTime() const noexcept
True if this is a write interval.
Definition: TimeStateI.H:73
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
label nFaces() const noexcept
Number of mesh faces.
defineTypeNameAndDebug(kinematicThinFilm, 0)
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
scalar definedPatchRadii_
List of radii for patches - if patch not defined, radius calculated based on mesh geometry...
Base class for film injection models, handling mass transfer from the film.
areaTensorField gradNHat_
Gradient of surface normals.
Macros for easy insertion into run-time selection tables.
tmp< scalarField > calcCosAngle(const edgeScalarField &phi) const
Calculate the cosine of the angle between gravity vector and cell out flow direction.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
tmp< areaScalarField > calcInvR1(const areaVectorField &U, const scalarField &calcCosAngle) const
Calculate local (inverse) radius of curvature.
dimensionedScalar pos(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
addToRunTimeSelectionTable(liquidFilmBase, kinematicThinFilm, dictionary)
const areaVectorField & Uf() const
Access const reference Uf.
const faMesh & regionMesh() const
Return the region mesh database.
const edgeScalarField & phi2s() const
Access continuity flux.
autoPtr< surfaceVectorField > Uf
const uniformDimensionedVectorField & g
const dimensionSet dimForce
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
int debug
Static debugging option.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
U
Definition: pEqn.H:72
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual const areaScalarField & rho() const =0
Access const reference rho.
Nothing to be read.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
An indexed form of CGAL::Triangulation_face_base_2<K> used to keep track of the vertices in the trian...
Definition: indexedFace.H:43
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:51
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:80
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity
const liquidFilmBase & film() const
Return const access to the film surface film model.