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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "curvatureSeparation.H"
30 #include "fvMesh.H"
31 #include "Time.H"
32 #include "volFields.H"
33 #include "kinematicSingleLayer.H"
34 #include "surfaceInterpolate.H"
35 #include "fvcDiv.H"
36 #include "fvcGrad.H"
37 #include "stringListOps.H"
38 #include "cyclicPolyPatch.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 namespace regionModels
45 {
46 namespace surfaceFilmModels
47 {
48 
49 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
50 
51 defineTypeNameAndDebug(curvatureSeparation, 0);
53 (
54  injectionModel,
55  curvatureSeparation,
56  dictionary
57 );
58 
59 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
60 
61 tmp<volScalarField> curvatureSeparation::calcInvR1
62 (
63  const volVectorField& U
64 ) const
65 {
66  // method 1
67 /*
68  tmp<volScalarField> tinvR1
69  (
70  new volScalarField("invR1", fvc::div(film().nHat()))
71  );
72 */
73 
74  // method 2
75  dimensionedScalar smallU("smallU", dimVelocity, ROOTVSMALL);
76  volVectorField UHat(U/(mag(U) + smallU));
77  tmp<volScalarField> tinvR1
78  (
79  new volScalarField("invR1", UHat & (UHat & gradNHat_))
80  );
81 
82 
83  scalarField& invR1 = tinvR1.ref().primitiveFieldRef();
84 
85  // apply defined patch radii
86  const scalar rMin = 1e-6;
87  const fvMesh& mesh = film().regionMesh();
88  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
89  forAll(definedPatchRadii_, i)
90  {
91  label patchi = definedPatchRadii_[i].first();
92  scalar definedInvR1 = 1.0/max(rMin, definedPatchRadii_[i].second());
93  UIndirectList<scalar>(invR1, pbm[patchi].faceCells()) = definedInvR1;
94  }
95 
96  // filter out large radii
97  const scalar rMax = 1e6;
98  forAll(invR1, i)
99  {
100  if (mag(invR1[i]) < 1/rMax)
101  {
102  invR1[i] = -1.0;
103  }
104  }
105 
106  if (debug && mesh.time().writeTime())
107  {
108  tinvR1().write();
109  }
110 
111  return tinvR1;
112 }
113 
114 
116 (
117  const surfaceScalarField& phi
118 ) const
119 {
120  const fvMesh& mesh = film().regionMesh();
121  const vectorField nf(mesh.Sf()/mesh.magSf());
122  const labelUList& own = mesh.owner();
123  const labelUList& nbr = mesh.neighbour();
124 
125  scalarField phiMax(mesh.nCells(), -GREAT);
126  scalarField cosAngle(mesh.nCells(), Zero);
127  forAll(nbr, facei)
128  {
129  label cellO = own[facei];
130  label cellN = nbr[facei];
131 
132  if (phi[facei] > phiMax[cellO])
133  {
134  phiMax[cellO] = phi[facei];
135  cosAngle[cellO] = -gHat_ & nf[facei];
136  }
137  if (-phi[facei] > phiMax[cellN])
138  {
139  phiMax[cellN] = -phi[facei];
140  cosAngle[cellN] = -gHat_ & -nf[facei];
141  }
142  }
143 
144  forAll(phi.boundaryField(), patchi)
145  {
146  const fvsPatchScalarField& phip = phi.boundaryField()[patchi];
147  const fvPatch& pp = phip.patch();
148  const labelList& faceCells = pp.faceCells();
149  const vectorField nf(pp.nf());
150  forAll(phip, i)
151  {
152  label celli = faceCells[i];
153  if (phip[i] > phiMax[celli])
154  {
155  phiMax[celli] = phip[i];
156  cosAngle[celli] = -gHat_ & nf[i];
157  }
158  }
159  }
160 /*
161  // correction for cyclics - use cyclic pairs' face normal instead of
162  // local face normal
163  const fvBoundaryMesh& pbm = mesh.boundary();
164  forAll(phi.boundaryField(), patchi)
165  {
166  if (isA<cyclicPolyPatch>(pbm[patchi]))
167  {
168  const scalarField& phip = phi.boundaryField()[patchi];
169  const vectorField nf(pbm[patchi].nf());
170  const labelList& faceCells = pbm[patchi].faceCells();
171  const label sizeBy2 = pbm[patchi].size()/2;
172 
173  for (label face0=0; face0<sizeBy2; face0++)
174  {
175  label face1 = face0 + sizeBy2;
176  label cell0 = faceCells[face0];
177  label cell1 = faceCells[face1];
178 
179  // flux leaving half 0, entering half 1
180  if (phip[face0] > phiMax[cell0])
181  {
182  phiMax[cell0] = phip[face0];
183  cosAngle[cell0] = -gHat_ & -nf[face1];
184  }
185 
186  // flux leaving half 1, entering half 0
187  if (-phip[face1] > phiMax[cell1])
188  {
189  phiMax[cell1] = -phip[face1];
190  cosAngle[cell1] = -gHat_ & nf[face0];
191  }
192  }
193  }
194  }
195 */
196  // checks
197  if (debug && mesh.time().writeTime())
198  {
199  volScalarField volCosAngle
200  (
201  IOobject
202  (
203  "cosAngle",
204  mesh.time().timeName(),
205  mesh,
207  ),
208  mesh,
210  zeroGradientFvPatchScalarField::typeName
211  );
212  volCosAngle.primitiveFieldRef() = cosAngle;
213  volCosAngle.correctBoundaryConditions();
214  volCosAngle.write();
215  }
216 
217  return max(min(cosAngle, scalar(1)), scalar(-1));
218 }
219 
220 
221 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
222 
223 curvatureSeparation::curvatureSeparation
224 (
226  const dictionary& dict
227 )
228 :
229  injectionModel(type(), film, dict),
230  gradNHat_(fvc::grad(film.nHat())),
231  deltaByR1Min_(coeffDict_.lookupOrDefault<scalar>("deltaByR1Min", 0)),
232  definedPatchRadii_(),
233  magG_(mag(film.g().value())),
234  gHat_(Zero)
235 {
236  if (magG_ < ROOTVSMALL)
237  {
239  << "Acceleration due to gravity must be non-zero"
240  << exit(FatalError);
241  }
242 
243  gHat_ = film.g().value()/magG_;
244 
245  List<Tuple2<word, scalar>> prIn(coeffDict_.lookup("definedPatchRadii"));
246  const wordList& allPatchNames = film.regionMesh().boundaryMesh().names();
247 
248  DynamicList<Tuple2<label, scalar>> prData(allPatchNames.size());
249 
250  labelHashSet uniquePatchIDs;
251 
252  forAllReverse(prIn, i)
253  {
254  labelList patchIDs = findIndices(allPatchNames, prIn[i].first());
255  forAll(patchIDs, j)
256  {
257  const label patchi = patchIDs[j];
258 
259  if (!uniquePatchIDs.found(patchi))
260  {
261  const scalar radius = prIn[i].second();
262  prData.append(Tuple2<label, scalar>(patchi, radius));
263 
264  uniquePatchIDs.insert(patchi);
265  }
266  }
267  }
268 
269  definedPatchRadii_.transfer(prData);
270 }
271 
272 
273 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
274 
276 {}
277 
278 
279 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
280 
282 (
283  scalarField& availableMass,
284  scalarField& massToInject,
285  scalarField& diameterToInject
286 )
287 {
288  const kinematicSingleLayer& film =
289  refCast<const kinematicSingleLayer>(this->film());
290  const fvMesh& mesh = film.regionMesh();
291 
292  const volScalarField& delta = film.delta();
293  const volVectorField& U = film.U();
294  const surfaceScalarField& phi = film.phi();
295  const volScalarField& rho = film.rho();
296  const scalarField magSqrU(magSqr(film.U()));
297  const volScalarField& sigma = film.sigma();
298 
299  const scalarField invR1(calcInvR1(U));
300  const scalarField cosAngle(calcCosAngle(phi));
301 
302  // calculate force balance
303  const scalar Fthreshold = 1e-10;
304  scalarField Fnet(mesh.nCells(), Zero);
305  scalarField separated(mesh.nCells(), Zero);
306  forAll(invR1, i)
307  {
308  if ((invR1[i] > 0) && (delta[i]*invR1[i] > deltaByR1Min_))
309  {
310  scalar R1 = 1.0/(invR1[i] + ROOTVSMALL);
311  scalar R2 = R1 + delta[i];
312 
313  // inertial force
314  scalar Fi = -delta[i]*rho[i]*magSqrU[i]*72.0/60.0*invR1[i];
315 
316  // body force
317  scalar Fb =
318  - 0.5*rho[i]*magG_*invR1[i]*(sqr(R1) - sqr(R2))*cosAngle[i];
319 
320  // surface force
321  scalar Fs = sigma[i]/R2;
322 
323  Fnet[i] = Fi + Fb + Fs;
324 
325  if (Fnet[i] + Fthreshold < 0)
326  {
327  separated[i] = 1.0;
328  }
329  }
330  }
331 
332  // inject all available mass
333  massToInject = separated*availableMass;
334  diameterToInject = separated*delta;
335  availableMass -= separated*availableMass;
336 
337  addToInjectedMass(sum(separated*availableMass));
338 
339  if (debug && mesh.time().writeTime())
340  {
341  volScalarField volFnet
342  (
343  IOobject
344  (
345  "Fnet",
346  mesh.time().timeName(),
347  mesh,
349  ),
350  mesh,
352  zeroGradientFvPatchScalarField::typeName
353  );
354  volFnet.primitiveFieldRef() = Fnet;
355  volFnet.correctBoundaryConditions();
356  volFnet.write();
357  }
358 
360 }
361 
362 
363 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
364 
365 } // End namespace surfaceFilmModels
366 } // End namespace regionModels
367 } // End namespace Foam
368 
369 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
volFields.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::regionModels::surfaceFilmModels::curvatureSeparation::~curvatureSeparation
virtual ~curvatureSeparation()
Destructor.
Definition: curvatureSeparation.C:275
Foam::regionModels::surfaceFilmModels::injectionModel::correct
void correct()
Correct.
Definition: injectionModel.C:81
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
cyclicPolyPatch.H
Foam::tmp< volScalarField >
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:57
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::U
virtual const volVectorField & U() const
Return the film velocity [m/s].
Definition: kinematicSingleLayer.C:981
Foam::regionModels::surfaceFilmModels::curvatureSeparation::calcInvR1
tmp< volScalarField > calcInvR1(const volVectorField &U) const
Calculate local (inverse) radius of curvature.
Definition: curvatureSeparation.C:62
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
fvcDiv.H
Calculate the divergence of the given field.
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer
Kinematic form of single-cell layer surface film model.
Definition: kinematicSingleLayer.H:66
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:182
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:404
Foam::dimForce
const dimensionSet dimForce
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Foam::HashSet< label, Hash< label > >
rho
rho
Definition: readInitialConditions.H:96
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyBoundaryMesh::names
wordList names() const
Return a list of patch names.
Definition: polyBoundaryMesh.C:593
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel::g
const dimensionedVector & g() const
Return the acceleration due to gravity.
Definition: surfaceFilmRegionModelI.H:41
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
CGAL::indexedFace
An indexed form of CGAL::Triangulation_face_base_2<K> used to keep track of the vertices in the trian...
Definition: indexedFace.H:50
Foam::fvPatch::nf
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:120
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::label
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:62
Foam::Field< scalar >
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel
Base class for surface film models.
Definition: surfaceFilmRegionModel.H:55
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:63
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::phi
virtual const surfaceScalarField & phi() const
Return the film flux [kg.m/s].
Definition: kinematicSingleLayer.C:1005
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::regionModels::surfaceFilmModels::injectionModel
Base class for film injection models, handling mass transfer from the film.
Definition: injectionModel.H:58
kinematicSingleLayer.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:735
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::sigma
const volScalarField & sigma() const
Return const access to the surface tension [kg/s2].
Definition: kinematicSingleLayerI.H:79
Time.H
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::fvsPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvsPatchField.H:281
Foam::List< label >
fvcGrad.H
Calculate the gradient of the given field.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::regionModels::singleLayerRegion::nHat
virtual const volVectorField & nHat() const
Return the patch normal vectors.
Definition: singleLayerRegion.C:209
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::delta
const volScalarField & delta() const
Return const access to the film thickness [m].
Definition: kinematicSingleLayerI.H:85
Foam::UList< label >
Foam::DynamicList::transfer
void transfer(List< T > &lst)
Transfer contents of the argument List into this.
Definition: DynamicListI.H:426
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:89
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:303
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:109
stringListOps.H
Operations on lists of strings.
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::Tuple2
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition: Tuple2.H:57
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rho
virtual const volScalarField & rho() const
Return the film density [kg/m3].
Definition: kinematicSingleLayer.C:1011
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
curvatureSeparation.H
Foam::regionModels::surfaceFilmModels::curvatureSeparation::calcCosAngle
tmp< scalarField > calcCosAngle(const surfaceScalarField &phi) const
Calculate the cosine of the angle between gravity vector and.
Definition: curvatureSeparation.C:116