curvatureSeparation.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "curvatureSeparation.H"
27 #include "fvcGrad.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace filmEjectionModels
36 {
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
42 (
46 );
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
50 tmp<scalarField> curvatureSeparation::calcInvR1
51 (
52  const volVectorField& U
53 ) const
54 {
55  const vectorField UHat(U.field()/(mag(U.field()) + rootVSmall));
56 
57  tmp<scalarField> tinvR1(-(UHat & (UHat & gradNHat_)));
58  scalarField& invR1 = tinvR1.ref();
59 
60  const scalar rMin = 1e-6;
61  const fvMesh& mesh = film_.mesh;
62  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
63 
64  forAll(patchRadii_, i)
65  {
66  const label patchi = patchRadii_[i].first();
67  const scalar definedInvR1 = 1/max(rMin, patchRadii_[i].second());
68  UIndirectList<scalar>(invR1, pbm[patchi].faceCells()) = definedInvR1;
69  }
70 
71  // Filter out large radii
72  const scalar rMax = 1e6;
73 
74  forAll(invR1, i)
75  {
76  if (mag(invR1[i]) < 1/rMax)
77  {
78  invR1[i] = -1;
79  }
80  }
81 
82  return tinvR1;
83 }
84 
85 
86 tmp<scalarField> curvatureSeparation::calcCosAngle
87 (
88  const surfaceScalarField& alphaRhoPhi
89 ) const
90 {
91  const fvMesh& mesh = film_.mesh;
92 
93  const scalar magg(mag(film_.g.value()));
94  const vector gHat(film_.g.value()/magg);
95 
96  const vectorField nf(mesh.Sf()/mesh.magSf());
97  const labelUList& own = mesh.owner();
98  const labelUList& nbr = mesh.neighbour();
99 
100  scalarField alphaRhoPhiMax(mesh.nCells(), -great);
101  tmp<scalarField> tcosAngle
102  (
103  new scalarField(mesh.nCells(), 0)
104  );
105  scalarField& cosAngle = tcosAngle.ref();
106 
107  forAll(nbr, facei)
108  {
109  const label cellO = own[facei];
110  const label cellN = nbr[facei];
111 
112  if (alphaRhoPhi[facei] > alphaRhoPhiMax[cellO])
113  {
114  alphaRhoPhiMax[cellO] = alphaRhoPhi[facei];
115  cosAngle[cellO] = -gHat & nf[facei];
116  }
117  if (-alphaRhoPhi[facei] > alphaRhoPhiMax[cellN])
118  {
119  alphaRhoPhiMax[cellN] = -alphaRhoPhi[facei];
120  cosAngle[cellN] = -gHat & -nf[facei];
121  }
122  }
123 
124  forAll(alphaRhoPhi.boundaryField(), patchi)
125  {
126  const fvsPatchScalarField& alphaRhoPhip =
127  alphaRhoPhi.boundaryField()[patchi];
128 
129  const fvPatch& pp = alphaRhoPhip.patch();
130 
131  if (pp.coupled())
132  {
133  const labelList& faceCells = pp.faceCells();
134  const vectorField nf(pp.nf());
135 
136  forAll(alphaRhoPhip, i)
137  {
138  const label celli = faceCells[i];
139 
140  if (alphaRhoPhip[i] > alphaRhoPhiMax[celli])
141  {
142  alphaRhoPhiMax[celli] = alphaRhoPhip[i];
143  cosAngle[celli] = -gHat & nf[i];
144  }
145  }
146  }
147  }
148 
149  return tcosAngle;
150 }
151 
152 
153 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
154 
156 (
157  const dictionary& dict,
158  const solvers::isothermalFilm& film
159 )
160 :
161  ejectionModel(dict, film),
162  gradNHat_(fvc::grad(film_.nHat)),
163  deltaByR1Min_
164  (
165  dict.optionalSubDict(typeName + "Coeffs")
166  .lookupOrDefault<scalar>("deltaByR1Min", scalar(0))
167  ),
168  deltaStable_
169  (
170  dict.optionalSubDict(typeName + "Coeffs")
171  .lookupOrDefault("deltaStable", scalar(0))
172  )
173 {
174  const List<Tuple2<word, scalar>> prIn
175  (
176  dict.lookupOrDefault("patchRadii", List<Tuple2<word, scalar>>::null())
177  );
178  const wordList& allPatchNames = film_.mesh.boundaryMesh().names();
179 
180  DynamicList<Tuple2<label, scalar>> prData(allPatchNames.size());
181 
182  labelHashSet uniquePatchIDs;
183 
184  forAllReverse(prIn, i)
185  {
186  labelList patchIDs = findStrings(prIn[i].first(), allPatchNames);
187  forAll(patchIDs, j)
188  {
189  const label patchi = patchIDs[j];
190 
191  if (!uniquePatchIDs.found(patchi))
192  {
193  const scalar radius = prIn[i].second();
194  prData.append(Tuple2<label, scalar>(patchi, radius));
195 
196  uniquePatchIDs.insert(patchi);
197  }
198  }
199  }
200 
201  patchRadii_.transfer(prData);
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
206 
208 {}
209 
210 
211 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
212 
214 {
215  const scalarField& delta = film_.delta();
216  const scalarField& rho = film_.rho();
217  const vectorField& U = film_.U();
218 
219  const tmp<volScalarField> tsigma = film_.sigma();
220  const volScalarField::Internal& sigma = tsigma();
221 
222  const scalarField invR1(calcInvR1(film_.U));
223  const scalarField cosAngle(calcCosAngle(film_.alphaRhoPhi));
224 
225  const scalar magg(mag(film_.g.value()));
226 
227  const scalar deltaT = film_.mesh.time().deltaTValue();
228 
229  // Calculate force balance
230  const scalar Fthreshold = 1e-10;
231 
232  forAll(delta, celli)
233  {
234  rate_[celli] = 0;
235  diameter_[celli] = 0;
236 
237  if
238  (
239  delta[celli] > deltaStable_
240  && invR1[celli] > 0
241  && delta[celli]*invR1[celli] > deltaByR1Min_
242  )
243  {
244  const scalar R1 = 1/(invR1[celli] + rootVSmall);
245  const scalar R2 = R1 + delta[celli];
246 
247  // Inertial force
248  const scalar Fi =
249  -delta[celli]*rho[celli]
250  *magSqr(U[celli])*72/60*invR1[celli];
251 
252  // Body force
253  const scalar Fb =
254  -0.5*rho[celli]*magg
255  *invR1[celli]*(sqr(R1) - sqr(R2))*cosAngle[celli];
256 
257  // Surface tension force
258  const scalar Fs = sigma[celli]/R2;
259 
260  if (Fi + Fb + Fs + Fthreshold < 0)
261  {
262  // Calculate rate of separation
263  rate_[celli] =
264  (delta[celli] - deltaStable_)/(deltaT*delta[celli]);
265 
266  // Set the separated droplet diameter to the film thickness
267  diameter_[celli] = delta[celli];
268  }
269  }
270  }
271 }
272 
273 
275 {
277  gradNHat_ = fvc::grad(film_.nHat);
278 }
279 
280 
282 {
284  gradNHat_ = fvc::grad(film_.nHat);
285 }
286 
287 
289 {
291  gradNHat_ = fvc::grad(film_.nHat);
292 }
293 
294 
296 {
298  gradNHat_ = fvc::grad(film_.nHat);
299 
300  return true;
301 }
302 
303 
304 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305 
306 } // End namespace filmEjectionModels
307 } // End namespace Foam
308 
309 // ************************************************************************* //
scalar delta
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
Generic GeometricField class.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:63
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
virtual bool movePoints()
Update for mesh motion.
Definition: ejectionModel.C:99
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: ejectionModel.C:75
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: ejectionModel.C:91
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: ejectionModel.C:83
Curvature induced separation film to cloud ejection transfer model.
virtual bool movePoints()
Update for mesh motion.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
curvatureSeparation(const dictionary &dict, const solvers::isothermalFilm &film)
Construct from dictionary and film model.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Foam::polyBoundaryMesh.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Solver module for flow of compressible isothermal liquid films.
A class for managing temporary objects.
Definition: tmp.H:55
Calculate the gradient of the given field.
label patchi
U
Definition: pEqn.H:72
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
addToRunTimeSelectionTable(ejectionModel, BrunDripping, dictionary)
defineTypeNameAndDebug(BrunDripping, 0)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:105
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
dimensionedSymmTensor sqr(const dimensionedVector &dv)
SurfaceField< scalar > surfaceScalarField
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
UList< label > labelUList
Definition: UList.H:65
dimensioned< scalar > magSqr(const dimensioned< Type > &)
fvsPatchField< scalar > fvsPatchScalarField
dictionary dict