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-2024 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
56  (
57  U.primitiveField()/(mag(U.primitiveField()) + rootVSmall)
58  );
59 
60  tmp<scalarField> tinvR1(-(UHat & (UHat & gradNHat_)));
61  scalarField& invR1 = tinvR1.ref();
62 
63  const scalar rMin = 1e-6;
64  const fvMesh& mesh = film_.mesh;
65  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
66 
67  forAll(patchRadii_, i)
68  {
69  const label patchi = patchRadii_[i].first();
70  const scalar definedInvR1 = 1/max(rMin, patchRadii_[i].second());
71  UIndirectList<scalar>(invR1, pbm[patchi].faceCells()) = definedInvR1;
72  }
73 
74  // Filter out large radii
75  const scalar rMax = 1e6;
76 
77  forAll(invR1, i)
78  {
79  if (mag(invR1[i]) < 1/rMax)
80  {
81  invR1[i] = -1;
82  }
83  }
84 
85  return tinvR1;
86 }
87 
88 
89 tmp<scalarField> curvatureSeparation::calcCosAngle
90 (
91  const surfaceScalarField& alphaRhoPhi
92 ) const
93 {
94  const fvMesh& mesh = film_.mesh;
95 
96  const scalar magg(mag(film_.g.value()));
97  const vector gHat(film_.g.value()/magg);
98 
99  const vectorField nf(mesh.Sf()/mesh.magSf());
100  const labelUList& own = mesh.owner();
101  const labelUList& nbr = mesh.neighbour();
102 
103  scalarField alphaRhoPhiMax(mesh.nCells(), -great);
104  tmp<scalarField> tcosAngle
105  (
106  new scalarField(mesh.nCells(), 0)
107  );
108  scalarField& cosAngle = tcosAngle.ref();
109 
110  forAll(nbr, facei)
111  {
112  const label cellO = own[facei];
113  const label cellN = nbr[facei];
114 
115  if (alphaRhoPhi[facei] > alphaRhoPhiMax[cellO])
116  {
117  alphaRhoPhiMax[cellO] = alphaRhoPhi[facei];
118  cosAngle[cellO] = -gHat & nf[facei];
119  }
120  if (-alphaRhoPhi[facei] > alphaRhoPhiMax[cellN])
121  {
122  alphaRhoPhiMax[cellN] = -alphaRhoPhi[facei];
123  cosAngle[cellN] = -gHat & -nf[facei];
124  }
125  }
126 
127  forAll(alphaRhoPhi.boundaryField(), patchi)
128  {
129  const fvsPatchScalarField& alphaRhoPhip =
130  alphaRhoPhi.boundaryField()[patchi];
131 
132  const fvPatch& pp = alphaRhoPhip.patch();
133 
134  if (pp.coupled())
135  {
136  const labelList& faceCells = pp.faceCells();
137  const vectorField nf(pp.nf());
138 
139  forAll(alphaRhoPhip, i)
140  {
141  const label celli = faceCells[i];
142 
143  if (alphaRhoPhip[i] > alphaRhoPhiMax[celli])
144  {
145  alphaRhoPhiMax[celli] = alphaRhoPhip[i];
146  cosAngle[celli] = -gHat & nf[i];
147  }
148  }
149  }
150  }
151 
152  return tcosAngle;
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
157 
159 (
160  const dictionary& dict,
161  const solvers::isothermalFilm& film
162 )
163 :
164  ejectionModel(dict, film),
165  gradNHat_(fvc::grad(film_.nHat)),
166  deltaByR1Min_
167  (
168  dict.optionalSubDict(typeName + "Coeffs")
169  .lookupOrDefault<scalar>("deltaByR1Min", scalar(0))
170  ),
171  deltaStable_
172  (
173  dict.optionalSubDict(typeName + "Coeffs")
174  .lookupOrDefault("deltaStable", scalar(0))
175  )
176 {
177  const List<Tuple2<word, scalar>> prIn
178  (
179  dict.lookupOrDefault("patchRadii", List<Tuple2<word, scalar>>::null())
180  );
181  const wordList& allPatchNames = film_.mesh.boundaryMesh().names();
182 
183  DynamicList<Tuple2<label, scalar>> prData(allPatchNames.size());
184 
185  labelHashSet uniquePatchIDs;
186 
187  forAllReverse(prIn, i)
188  {
189  labelList patchIDs = findStrings(prIn[i].first(), allPatchNames);
190  forAll(patchIDs, j)
191  {
192  const label patchi = patchIDs[j];
193 
194  if (!uniquePatchIDs.found(patchi))
195  {
196  const scalar radius = prIn[i].second();
197  prData.append(Tuple2<label, scalar>(patchi, radius));
198 
199  uniquePatchIDs.insert(patchi);
200  }
201  }
202  }
203 
204  patchRadii_.transfer(prData);
205 }
206 
207 
208 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
209 
211 {}
212 
213 
214 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
215 
217 {
218  const scalarField& delta = film_.delta();
219  const scalarField& rho = film_.rho();
220  const vectorField& U = film_.U();
221 
222  const tmp<volScalarField> tsigma = film_.sigma();
223  const volScalarField::Internal& sigma = tsigma();
224 
225  const scalarField invR1(calcInvR1(film_.U));
226  const scalarField cosAngle(calcCosAngle(film_.alphaRhoPhi));
227 
228  const scalar magg(mag(film_.g.value()));
229 
230  const scalar deltaT = film_.mesh.time().deltaTValue();
231 
232  // Calculate force balance
233  const scalar Fthreshold = 1e-10;
234 
235  forAll(delta, celli)
236  {
237  rate_[celli] = 0;
238  diameter_[celli] = 0;
239 
240  if
241  (
242  delta[celli] > deltaStable_
243  && invR1[celli] > 0
244  && delta[celli]*invR1[celli] > deltaByR1Min_
245  )
246  {
247  const scalar R1 = 1/(invR1[celli] + rootVSmall);
248  const scalar R2 = R1 + delta[celli];
249 
250  // Inertial force
251  const scalar Fi =
252  -delta[celli]*rho[celli]
253  *magSqr(U[celli])*72/60*invR1[celli];
254 
255  // Body force
256  const scalar Fb =
257  -0.5*rho[celli]*magg
258  *invR1[celli]*(sqr(R1) - sqr(R2))*cosAngle[celli];
259 
260  // Surface tension force
261  const scalar Fs = sigma[celli]/R2;
262 
263  if (Fi + Fb + Fs + Fthreshold < 0)
264  {
265  // Calculate rate of separation
266  rate_[celli] =
267  (delta[celli] - deltaStable_)/(deltaT*delta[celli]);
268 
269  // Set the separated droplet diameter to the film thickness
270  diameter_[celli] = delta[celli];
271  }
272  }
273  }
274 }
275 
276 
278 {
280  gradNHat_ = fvc::grad(film_.nHat);
281 }
282 
283 
285 {
287  gradNHat_ = fvc::grad(film_.nHat);
288 }
289 
290 
292 {
294  gradNHat_ = fvc::grad(film_.nHat);
295 }
296 
297 
299 {
301  gradNHat_ = fvc::grad(film_.nHat);
302 
303  return true;
304 }
305 
306 
307 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308 
309 } // End namespace filmEjectionModels
310 } // End namespace Foam
311 
312 // ************************************************************************* //
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:109
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
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:66
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
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:99
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:441
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:404
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:106
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