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 "fvMesh.H"
32 #include "Time.H"
33 #include "volFields.H"
34 #include "kinematicSingleLayer.H"
35 #include "surfaceInterpolate.H"
36 #include "fvcDiv.H"
37 #include "fvcGrad.H"
38 #include "stringListOps.H"
39 #include "cyclicPolyPatch.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace regionModels
46 {
47 namespace surfaceFilmModels
48 {
49 
50 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
51 
52 defineTypeNameAndDebug(curvatureSeparation, 0);
54 (
55  injectionModel,
58 );
59 
60 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
61 
63 (
64  const volVectorField& U
65 ) const
66 {
67  // method 1
68 /*
69  tmp<volScalarField> tinvR1
70  (
71  new volScalarField("invR1", fvc::div(film().nHat()))
72  );
73 */
74 
75  // method 2
76  dimensionedScalar smallU("smallU", dimVelocity, ROOTVSMALL);
77  volVectorField UHat(U/(mag(U) + smallU));
78  tmp<volScalarField> tinvR1
79  (
80  new volScalarField("invR1", UHat & (UHat & gradNHat_))
81  );
82 
83 
84  scalarField& invR1 = tinvR1.ref().primitiveFieldRef();
85 
86  // apply defined patch radii
87  const scalar rMin = 1e-6;
88  const fvMesh& mesh = film().regionMesh();
89  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
91  {
92  label patchi = definedPatchRadii_[i].first();
93  scalar definedInvR1 = 1.0/max(rMin, definedPatchRadii_[i].second());
94  UIndirectList<scalar>(invR1, pbm[patchi].faceCells()) = definedInvR1;
95  }
96 
97  // filter out large radii
98  const scalar rMax = 1e6;
99  forAll(invR1, i)
100  {
101  if (mag(invR1[i]) < 1/rMax)
102  {
103  invR1[i] = -1.0;
104  }
105  }
106 
107  if (debug && mesh.time().writeTime())
108  {
109  tinvR1().write();
110  }
111 
112  return tinvR1;
113 }
114 
115 
117 (
118  const surfaceScalarField& phi
119 ) const
120 {
121  const fvMesh& mesh = film().regionMesh();
122  const vectorField nf(mesh.Sf()/mesh.magSf());
123  const labelUList& own = mesh.owner();
124  const labelUList& nbr = mesh.neighbour();
125 
126  scalarField phiMax(mesh.nCells(), -GREAT);
127  scalarField cosAngle(mesh.nCells(), Zero);
128  forAll(nbr, facei)
129  {
130  label cellO = own[facei];
131  label cellN = nbr[facei];
132 
133  if (phi[facei] > phiMax[cellO])
134  {
135  phiMax[cellO] = phi[facei];
136  cosAngle[cellO] = -gHat_ & nf[facei];
137  }
138  if (-phi[facei] > phiMax[cellN])
139  {
140  phiMax[cellN] = -phi[facei];
141  cosAngle[cellN] = -gHat_ & -nf[facei];
142  }
143  }
144 
145  forAll(phi.boundaryField(), patchi)
146  {
147  const fvsPatchScalarField& phip = phi.boundaryField()[patchi];
148  const fvPatch& pp = phip.patch();
149  const labelList& faceCells = pp.faceCells();
150  const vectorField nf(pp.nf());
151  forAll(phip, i)
152  {
153  label celli = faceCells[i];
154  if (phip[i] > phiMax[celli])
155  {
156  phiMax[celli] = phip[i];
157  cosAngle[celli] = -gHat_ & nf[i];
158  }
159  }
160  }
161 /*
162  // correction for cyclics - use cyclic pairs' face normal instead of
163  // local face normal
164  const fvBoundaryMesh& pbm = mesh.boundary();
165  forAll(phi.boundaryField(), patchi)
166  {
167  if (isA<cyclicPolyPatch>(pbm[patchi]))
168  {
169  const scalarField& phip = phi.boundaryField()[patchi];
170  const vectorField nf(pbm[patchi].nf());
171  const labelList& faceCells = pbm[patchi].faceCells();
172  const label sizeBy2 = pbm[patchi].size()/2;
173 
174  for (label face0=0; face0<sizeBy2; face0++)
175  {
176  label face1 = face0 + sizeBy2;
177  label cell0 = faceCells[face0];
178  label cell1 = faceCells[face1];
179 
180  // flux leaving half 0, entering half 1
181  if (phip[face0] > phiMax[cell0])
182  {
183  phiMax[cell0] = phip[face0];
184  cosAngle[cell0] = -gHat_ & -nf[face1];
185  }
186 
187  // flux leaving half 1, entering half 0
188  if (-phip[face1] > phiMax[cell1])
189  {
190  phiMax[cell1] = -phip[face1];
191  cosAngle[cell1] = -gHat_ & nf[face0];
192  }
193  }
194  }
195  }
196 */
197  // checks
198  if (debug && mesh.time().writeTime())
199  {
200  volScalarField volCosAngle
201  (
202  IOobject
203  (
204  "cosAngle",
205  mesh.time().timeName(),
206  mesh,
208  ),
209  mesh,
211  zeroGradientFvPatchScalarField::typeName
212  );
213  volCosAngle.primitiveFieldRef() = cosAngle;
214  volCosAngle.correctBoundaryConditions();
215  volCosAngle.write();
216  }
217 
218  return max(min(cosAngle, scalar(1)), scalar(-1));
219 }
220 
221 
222 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
223 
224 curvatureSeparation::curvatureSeparation
225 (
227  const dictionary& dict
228 )
229 :
230  injectionModel(type(), film, dict),
231  gradNHat_(fvc::grad(film.nHat())),
232  deltaByR1Min_(coeffDict_.getOrDefault<scalar>("deltaByR1Min", 0)),
233  definedPatchRadii_(),
234  magG_(mag(film.g().value())),
235  gHat_(Zero)
236 {
237  if (magG_ < ROOTVSMALL)
238  {
240  << "Acceleration due to gravity must be non-zero"
241  << exit(FatalError);
242  }
243 
244  gHat_ = film.g().value()/magG_;
245 
246  List<Tuple2<word, scalar>> prIn(coeffDict_.lookup("definedPatchRadii"));
247  const wordList& allPatchNames = film.regionMesh().boundaryMesh().names();
248 
249  DynamicList<Tuple2<label, scalar>> prData(allPatchNames.size());
250 
251  labelHashSet uniquePatchIDs;
252 
253  forAllReverse(prIn, i)
254  {
255  labelList patchIDs = findIndices(allPatchNames, prIn[i].first());
256  forAll(patchIDs, j)
257  {
258  const label patchi = patchIDs[j];
259 
260  if (!uniquePatchIDs.found(patchi))
261  {
262  const scalar radius = prIn[i].second();
263  prData.append(Tuple2<label, scalar>(patchi, radius));
264 
265  uniquePatchIDs.insert(patchi);
266  }
267  }
268  }
270  definedPatchRadii_.transfer(prData);
271 }
272 
273 
274 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
275 
277 {}
278 
279 
280 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
281 
283 (
284  scalarField& availableMass,
285  scalarField& massToInject,
286  scalarField& diameterToInject
287 )
288 {
289  const kinematicSingleLayer& film =
290  refCast<const kinematicSingleLayer>(this->film());
291  const fvMesh& mesh = film.regionMesh();
292 
293  const volScalarField& delta = film.delta();
294  const volVectorField& U = film.U();
295  const surfaceScalarField& phi = film.phi();
296  const volScalarField& rho = film.rho();
297  const scalarField magSqrU(magSqr(film.U()));
298  const volScalarField& sigma = film.sigma();
299 
300  const scalarField invR1(calcInvR1(U));
301  const scalarField cosAngle(calcCosAngle(phi));
302 
303  // calculate force balance
304  const scalar Fthreshold = 1e-10;
305  scalarField Fnet(mesh.nCells(), Zero);
306  scalarField separated(mesh.nCells(), Zero);
307  forAll(invR1, i)
308  {
309  if ((invR1[i] > 0) && (delta[i]*invR1[i] > deltaByR1Min_))
310  {
311  scalar R1 = 1.0/(invR1[i] + ROOTVSMALL);
312  scalar R2 = R1 + delta[i];
313 
314  // inertial force
315  scalar Fi = -delta[i]*rho[i]*magSqrU[i]*72.0/60.0*invR1[i];
316 
317  // body force
318  scalar Fb =
319  - 0.5*rho[i]*magG_*invR1[i]*(sqr(R1) - sqr(R2))*cosAngle[i];
320 
321  // surface force
322  scalar Fs = sigma[i]/R2;
323 
324  Fnet[i] = Fi + Fb + Fs;
325 
326  if (Fnet[i] + Fthreshold < 0)
327  {
328  separated[i] = 1.0;
329  }
330  }
331  }
332 
333  // inject all available mass
334  massToInject = separated*availableMass;
335  diameterToInject = separated*delta;
336  availableMass -= separated*availableMass;
337 
338  addToInjectedMass(sum(separated*availableMass));
339 
340  if (debug && mesh.time().writeTime())
341  {
342  volScalarField volFnet
343  (
344  IOobject
345  (
346  "Fnet",
347  mesh.time().timeName(),
348  mesh,
350  ),
351  mesh,
353  zeroGradientFvPatchScalarField::typeName
354  );
355  volFnet.primitiveFieldRef() = Fnet;
356  volFnet.correctBoundaryConditions();
357  volFnet.write();
358  }
359 
361 }
362 
363 
364 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
365 
366 } // End namespace surfaceFilmModels
367 } // End namespace regionModels
368 } // End namespace Foam
369 
370 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
List< Tuple2< label, scalar > > definedPatchRadii_
List of radii for patches - if patch not defined, radius.
scalar delta
fvsPatchField< scalar > fvsPatchScalarField
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream. ...
Definition: dictionary.C:379
const Type & value() const noexcept
Return const reference to value.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
tmp< scalarField > calcCosAngle(const surfaceScalarField &phi) const
Calculate the cosine of the angle between gravity vector and.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Kinematic form of single-cell layer surface film model.
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:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: stringOps.H:55
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:491
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
bool found(const T &val, label pos=0) const
True if the value if found in the list.
Definition: UListI.H:258
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
Operations on lists of strings.
const dimensionSet dimless
Dimensionless.
tmp< volScalarField > calcInvR1(const volVectorField &U) const
Calculate local (inverse) radius of curvature.
const dictionary coeffDict_
Coefficients dictionary.
Definition: subModelBase.H:92
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:752
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
Calculate the gradient of the given field.
wordList names() const
Return a list of patch names.
virtual const volScalarField & delta() const =0
Return the film thickness [m].
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:26
virtual const volVectorField & U() const =0
Return the film velocity [m/s].
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
Calculate the divergence of the given field.
void addToInjectedMass(const scalar dMass)
Add to injected mass.
const dimensionedVector & g() const
Return the acceleration due to gravity.
const uniformDimensionedVectorField & g
const dimensionSet dimForce
int debug
Static debugging option.
Base class for film injection models, handling mass transfer from the film.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:57
Nothing to be read.
Field< vector > vectorField
Specialisation of Field<T> for vector.
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:429
List< label > labelList
A List of labels.
Definition: List.H:62
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
volTensorField gradNHat_
Gradient of surface normals.
virtual const volScalarField & sigma() const =0
Return the film surface tension [N/m].
defineTypeNameAndDebug(kinematicSingleLayer, 0)
scalar deltaByR1Min_
Minimum gravity driven film thickness (non-dimensionalised delta/R1)
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:49
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157
const dimensionSet dimVelocity