Loading...
Searching...
No Matches
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-2023 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 "curvatureSeparation.H"
31#include "fvMesh.H"
32#include "Time.H"
33#include "volFields.H"
35#include "surfaceInterpolate.H"
36#include "fvcDiv.H"
37#include "fvcGrad.H"
38#include "cyclicPolyPatch.H"
39
40// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41
42namespace Foam
43{
44namespace regionModels
46namespace surfaceFilmModels
47{
48
49// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
50
53(
54 injectionModel,
55 curvatureSeparation,
57);
58
59// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
60
62(
63 const volVectorField& U
64) const
65{
66 // method 1
67/*
68 auto tinvR1 = volScalarField::New("invR1", fvc::div(film().nHat()));
69*/
70
71 // method 2
72 dimensionedScalar smallU("smallU", dimVelocity, ROOTVSMALL);
73 volVectorField UHat(U/(mag(U) + smallU));
74
75 auto tinvR1 = volScalarField::New
76 (
77 "invR1",
79 UHat & (UHat & gradNHat_)
80 );
81
82 scalarField& invR1 = tinvR1.ref().primitiveFieldRef();
83
84 // apply defined patch radii
85 const scalar rMin = 1e-6;
86 const fvMesh& mesh = film().regionMesh();
87 const polyBoundaryMesh& pbm = mesh.boundaryMesh();
89 {
90 label patchi = definedPatchRadii_[i].first();
91 scalar definedInvR1 = 1.0/max(rMin, definedPatchRadii_[i].second());
92 UIndirectList<scalar>(invR1, pbm[patchi].faceCells()) = definedInvR1;
93 }
94
95 // filter out large radii
96 const scalar rMax = 1e6;
97 forAll(invR1, i)
98 {
99 if (mag(invR1[i]) < 1/rMax)
100 {
101 invR1[i] = -1.0;
102 }
103 }
104
105 if (debug && mesh.time().writeTime())
106 {
107 tinvR1().write();
108 }
109
110 return tinvR1;
111}
112
113
115(
117) const
118{
119 const fvMesh& mesh = film().regionMesh();
120 const vectorField nf(mesh.Sf()/mesh.magSf());
121 const labelUList& own = mesh.owner();
122 const labelUList& nbr = mesh.neighbour();
123
124 scalarField phiMax(mesh.nCells(), -GREAT);
125 scalarField cosAngle(mesh.nCells(), Zero);
126 forAll(nbr, facei)
127 {
128 label cellO = own[facei];
129 label cellN = nbr[facei];
130
131 if (phi[facei] > phiMax[cellO])
132 {
133 phiMax[cellO] = phi[facei];
134 cosAngle[cellO] = -gHat_ & nf[facei];
135 }
136 if (-phi[facei] > phiMax[cellN])
137 {
138 phiMax[cellN] = -phi[facei];
139 cosAngle[cellN] = -gHat_ & -nf[facei];
140 }
141 }
142
143 forAll(phi.boundaryField(), patchi)
144 {
145 const fvsPatchScalarField& phip = phi.boundaryField()[patchi];
146 const fvPatch& pp = phip.patch();
147 const labelUList& faceCells = pp.faceCells();
148 const vectorField nf(pp.nf());
149 forAll(phip, i)
150 {
151 label celli = faceCells[i];
152 if (phip[i] > phiMax[celli])
153 {
154 phiMax[celli] = phip[i];
155 cosAngle[celli] = -gHat_ & nf[i];
156 }
157 }
158 }
159/*
160 // correction for cyclics - use cyclic pairs' face normal instead of
161 // local face normal
162 const fvBoundaryMesh& pbm = mesh.boundary();
163 forAll(phi.boundaryField(), patchi)
164 {
165 if (isA<cyclicPolyPatch>(pbm[patchi]))
166 {
167 const scalarField& phip = phi.boundaryField()[patchi];
168 const vectorField nf(pbm[patchi].nf());
169 const labelUList& faceCells = pbm[patchi].faceCells();
170 const label sizeBy2 = pbm[patchi].size()/2;
171
172 for (label face0=0; face0<sizeBy2; face0++)
173 {
174 label face1 = face0 + sizeBy2;
175 label cell0 = faceCells[face0];
176 label cell1 = faceCells[face1];
177
178 // flux leaving half 0, entering half 1
179 if (phip[face0] > phiMax[cell0])
180 {
181 phiMax[cell0] = phip[face0];
182 cosAngle[cell0] = -gHat_ & -nf[face1];
183 }
184
185 // flux leaving half 1, entering half 0
186 if (-phip[face1] > phiMax[cell1])
187 {
188 phiMax[cell1] = -phip[face1];
189 cosAngle[cell1] = -gHat_ & nf[face0];
190 }
191 }
192 }
193 }
194*/
195 // checks
196 if (debug && mesh.time().writeTime())
197 {
198 auto tvolCosAngle = volScalarField::New
199 (
200 "cosAngle",
202 mesh,
205 );
206 auto& volCosAngle = tvolCosAngle.ref();
207
208 volCosAngle.primitiveFieldRef() = cosAngle;
209 volCosAngle.correctBoundaryConditions();
210 volCosAngle.write();
211 }
213 return clamp(cosAngle, scalarMinMax(-1, 1));
214}
215
216
217// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
218
219curvatureSeparation::curvatureSeparation
220(
222 const dictionary& dict
223)
224:
225 injectionModel(type(), film, dict),
226 gradNHat_(fvc::grad(film.nHat())),
227 deltaByR1Min_(coeffDict_.getOrDefault<scalar>("deltaByR1Min", 0)),
228 definedPatchRadii_(),
229 magG_(mag(film.g().value())),
230 gHat_(Zero)
231{
232 if (magG_ < ROOTVSMALL)
233 {
235 << "Acceleration due to gravity must be non-zero"
236 << exit(FatalError);
237 }
238
239 gHat_ = film.g().value()/magG_;
240
241 List<Tuple2<word, scalar>> prIn(coeffDict_.lookup("definedPatchRadii"));
242 const wordList& allPatchNames = film.regionMesh().boundaryMesh().names();
243
244 DynamicList<Tuple2<label, scalar>> prData(allPatchNames.size());
245
246 labelHashSet uniquePatchIDs;
247
248 forAllReverse(prIn, i)
249 {
250 labelList patchIDs = findIndices(allPatchNames, prIn[i].first());
251 forAll(patchIDs, j)
252 {
253 const label patchi = patchIDs[j];
254
255 if (!uniquePatchIDs.found(patchi))
256 {
257 const scalar radius = prIn[i].second();
258 prData.append(Tuple2<label, scalar>(patchi, radius));
259
260 uniquePatchIDs.insert(patchi);
261 }
262 }
263 }
265 definedPatchRadii_.transfer(prData);
266}
267
268
269// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
272{}
273
274
275// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
276
278(
279 scalarField& availableMass,
280 scalarField& massToInject,
281 scalarField& diameterToInject
282)
283{
286 const fvMesh& mesh = film.regionMesh();
287
288 const volScalarField& delta = film.delta();
289 const volVectorField& U = film.U();
290 const surfaceScalarField& phi = film.phi();
291 const volScalarField& rho = film.rho();
292 const scalarField magSqrU(magSqr(film.U()));
293 const volScalarField& sigma = film.sigma();
294
295 const scalarField invR1(calcInvR1(U));
296 const scalarField cosAngle(calcCosAngle(phi));
297
298 // calculate force balance
299 const scalar Fthreshold = 1e-10;
300 scalarField Fnet(mesh.nCells(), Zero);
301 scalarField separated(mesh.nCells(), Zero);
302 forAll(invR1, i)
303 {
304 if ((invR1[i] > 0) && (delta[i]*invR1[i] > deltaByR1Min_))
305 {
306 scalar R1 = 1.0/(invR1[i] + ROOTVSMALL);
307 scalar R2 = R1 + delta[i];
308
309 // inertial force
310 scalar Fi = -delta[i]*rho[i]*magSqrU[i]*72.0/60.0*invR1[i];
311
312 // body force
313 scalar Fb =
314 - 0.5*rho[i]*magG_*invR1[i]*(sqr(R1) - sqr(R2))*cosAngle[i];
315
316 // surface force
317 scalar Fs = sigma[i]/R2;
318
319 Fnet[i] = Fi + Fb + Fs;
320
321 if (Fnet[i] + Fthreshold < 0)
322 {
323 separated[i] = 1.0;
324 }
325 }
326 }
327
328 // inject all available mass
329 massToInject = separated*availableMass;
330 diameterToInject = separated*delta;
331 availableMass -= separated*availableMass;
332
333 addToInjectedMass(sum(separated*availableMass));
334
335 if (debug && mesh.time().writeTime())
336 {
337 auto tvolFnet = volScalarField::New
338 (
339 "Fnet",
341 mesh,
344 );
345 auto& volFnet = tvolFnet.ref();
346
347 volFnet.primitiveFieldRef() = Fnet;
348 volFnet.correctBoundaryConditions();
349 volFnet.write();
350 }
351
353}
354
355
356// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
357
358} // End namespace surfaceFilmModels
359} // End namespace regionModels
360} // End namespace Foam
361
362// ************************************************************************* //
CGAL::indexedFace< K > Fb
scalar delta
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
const polyBoundaryMesh & pbm
const uniformDimensionedVectorField & g
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
@ NO_REGISTER
Do not request registration (bool: false).
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition Tuple2.H:51
A List with indirect addressing. Like IndirectList but does not store addressing.
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
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:367
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
const fvPatch & patch() const noexcept
Return the patch.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
wordList names() const
Return a list of patch names.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const fvMesh & regionMesh() const
Return the region mesh database.
scalar deltaByR1Min_
Minimum gravity driven film thickness (non-dimensionalised delta/R1).
volTensorField gradNHat_
Gradient of surface normals.
tmp< volScalarField > calcInvR1(const volVectorField &U) const
Calculate local (inverse) radius of curvature.
tmp< scalarField > calcCosAngle(const surfaceScalarField &phi) const
Calculate the cosine of the angle between gravity vector and.
List< Tuple2< label, scalar > > definedPatchRadii_
List of radii for patches - if patch not defined, radius.
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
Base class for film injection models, handling mass transfer from the film.
void addToInjectedMass(const scalar dMass)
Add to injected mass.
Kinematic form of single-cell layer surface film model.
virtual const volScalarField & sigma() const =0
Return the film surface tension [N/m].
virtual const volVectorField & U() const =0
Return the film velocity [m/s].
virtual const volScalarField & rho() const =0
Return the film density [kg/m3].
const dimensionedVector & g() const
Return the acceleration due to gravity.
virtual const volScalarField & delta() const =0
Return the film thickness [m].
const dictionary coeffDict_
Coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace of functions to calculate explicit derivatives.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionSet dimForce
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvsPatchField< scalar > fvsPatchScalarField
dictionary dict
volScalarField & e
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition stdFoam.H:315