plicRDF.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) 2019-2020 DLR
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 "plicRDF.H"
29 #include "interpolationCellPoint.H"
30 #include "fvc.H"
31 #include "leastSquareGrad.H"
33 #include "profiling.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace reconstruction
40 {
41  defineTypeNameAndDebug(plicRDF, 0);
42  addToRunTimeSelectionTable(reconstructionSchemes,plicRDF, components);
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::reconstruction::plicRDF::interpolateNormal()
50 {
51  addProfilingInFunction(geometricVoF);
52  scalar dt = mesh_.time().deltaTValue();
53  zoneDistribute& exchangeFields = zoneDistribute::New(mesh_);
54 
55  leastSquareGrad<scalar> lsGrad("polyDegree1",mesh_.geometricD());
56 
57  exchangeFields.setUpCommforZone(interfaceCell_,false);
58 
59  Map<vector> mapCentre
60  (
61  exchangeFields.getDatafromOtherProc(interfaceCell_, centre_)
62  );
63  Map<vector> mapNormal
64  (
65  exchangeFields.getDatafromOtherProc(interfaceCell_, normal_)
66  );
67 
68  Map<vector> mapCC
69  (
70  exchangeFields.getDatafromOtherProc(interfaceCell_, mesh_.C())
71  );
72  Map<scalar> mapAlpha
73  (
74  exchangeFields.getDatafromOtherProc(interfaceCell_, alpha1_)
75  );
76 
77  DynamicField<vector > cellCentre(100);
78  DynamicField<scalar > alphaValues(100);
79 
80  DynamicList<vector> foundNormals(30);
81 
82  const labelListList& stencil = exchangeFields.getStencil();
83 
85  {
86  const label celli = interfaceLabels_[i];
87  vector estimatedNormal{Zero};
88  scalar weight{0};
89  foundNormals.clear();
90  for (const label gblIdx : stencil[celli])
91  {
92  vector n =
93  exchangeFields.getValue(normal_, mapNormal, gblIdx);
94  point p = mesh_.C()[celli]-U_[celli]*dt;
95  if (mag(n) != 0)
96  {
97  n /= mag(n);
98  vector centre =
99  exchangeFields.getValue(centre_, mapCentre, gblIdx);
100  vector distanceToIntSeg = (tensor::I- n*n) & (p - centre);
101  estimatedNormal += n /max(mag(distanceToIntSeg), SMALL);
102  weight += 1/max(mag(distanceToIntSeg), SMALL);
103  foundNormals.append(n);
104  }
105  }
106 
107  if (weight != 0 && mag(estimatedNormal) != 0)
108  {
109  estimatedNormal /= weight;
110  estimatedNormal /= mag(estimatedNormal);
111  }
112 
113  bool tooCoarse = false;
114 
115  if (foundNormals.size() > 1 && mag(estimatedNormal) != 0)
116  {
117  forAll(foundNormals, i)
118  {
119  // all have the length of 1
120  // to coarse if normal angle is bigger than 10 deg
121  if ((estimatedNormal & foundNormals[i]) <= 0.98)
122  {
123  tooCoarse = true;
124  }
125  }
126  }
127  else
128  {
129  tooCoarse = true;
130  }
131 
132  // if a normal was found and the interface is fine enough
133  // smallDist is always smallDist
134  if (mag(estimatedNormal) != 0 && !tooCoarse)
135  {
136  interfaceNormal_[i] = estimatedNormal;
137  }
138  else
139  {
140  cellCentre.clear();
141  alphaValues.clear();
142 
143  forAll(stencil[celli],i)
144  {
145  const label gblIdx = stencil[celli][i];
146  cellCentre.append
147  (
148  exchangeFields.getValue(mesh_.C(), mapCC, gblIdx)
149  );
150  alphaValues.append
151  (
152  exchangeFields.getValue(alpha1_, mapAlpha, gblIdx)
153  );
154  }
155  cellCentre -= mesh_.C()[celli];
156  interfaceNormal_[i] = lsGrad.grad(cellCentre, alphaValues);
157  }
158 
159  }
160 }
161 
162 void Foam::reconstruction::plicRDF::gradSurf(const volScalarField& phi)
163 {
164  addProfilingInFunction(geometricVoF);
165  leastSquareGrad<scalar> lsGrad("polyDegree1", mesh_.geometricD());
166  zoneDistribute& exchangeFields = zoneDistribute::New(mesh_);
167 
168  exchangeFields.setUpCommforZone(interfaceCell_, false);
169 
170  Map<vector> mapCC
171  (
172  exchangeFields.getDatafromOtherProc(interfaceCell_, mesh_.C())
173  );
174  Map<scalar> mapPhi
175  (
176  exchangeFields.getDatafromOtherProc(interfaceCell_, phi)
177  );
178 
179  DynamicField<vector> cellCentre(100);
180  DynamicField<scalar> phiValues(100);
181 
182  const labelListList& stencil = exchangeFields.getStencil();
183 
184  forAll(interfaceLabels_, i)
185  {
186  const label celli = interfaceLabels_[i];
187 
188  cellCentre.clear();
189  phiValues.clear();
190 
191  for (const label gblIdx : stencil[celli])
192  {
193  cellCentre.append
194  (
195  exchangeFields.getValue(mesh_.C(), mapCC, gblIdx)
196  );
197  phiValues.append
198  (
199  exchangeFields.getValue(phi, mapPhi, gblIdx)
200  );
201  }
202 
203  cellCentre -= mesh_.C()[celli];
204  interfaceNormal_[i] = lsGrad.grad(cellCentre, phiValues);
205  }
206 }
207 
208 
209 void Foam::reconstruction::plicRDF::setInitNormals(bool interpolate)
210 {
211  addProfilingInFunction(geometricVoF);
212  zoneDistribute& exchangeFields = zoneDistribute::New(mesh_);
213 
214  interfaceLabels_.clear();
215 
216  forAll(alpha1_, celli)
217  {
218  if (sIterPLIC_.isASurfaceCell(alpha1_[celli]))
219  {
220  interfaceCell_[celli] = true; // is set to false earlier
221  interfaceLabels_.append(celli);
222  }
223  }
224  interfaceNormal_.setSize(interfaceLabels_.size());
225 
226  RDF_.markCellsNearSurf(interfaceCell_, 1);
227  const boolList& nextToInterface_ = RDF_.nextToInterface();
228  exchangeFields.updateStencil(nextToInterface_);
229 
230  if (interpolate)
231  {
232  interpolateNormal();
233  }
234  else
235  {
236  gradSurf(alpha1_);
237  }
238 }
239 
240 
241 void Foam::reconstruction::plicRDF::calcResidual
242 (
243  List<normalRes>& normalResidual
244 )
245 {
246  addProfilingInFunction(geometricVoF);
247  zoneDistribute& exchangeFields = zoneDistribute::New(mesh_);
248  exchangeFields.setUpCommforZone(interfaceCell_,false);
249 
250  Map<vector> mapNormal
251  (
252  exchangeFields.getDatafromOtherProc(interfaceCell_, normal_)
253  );
254 
255  const labelListList& stencil = exchangeFields.getStencil();
256 
257  forAll(interfaceLabels_, i)
258  {
259  const label celli = interfaceLabels_[i];
260  if (mag(normal_[celli]) == 0 || mag(interfaceNormal_[i]) == 0)
261  {
262  normalResidual[i].celli = celli;
263  normalResidual[i].normalResidual = 0;
264  normalResidual[i].avgAngle = 0;
265  continue;
266  }
267 
268  scalar avgDiffNormal = 0;
269  scalar maxDiffNormal = GREAT;
270  scalar weight= 0;
271  const vector cellNormal = normal_[celli]/mag(normal_[celli]);
272  forAll(stencil[celli],j)
273  {
274  const label gblIdx = stencil[celli][j];
275  vector normal =
276  exchangeFields.getValue(normal_, mapNormal, gblIdx);
277 
278  if (mag(normal) != 0 && j != 0)
279  {
280  vector n = normal/mag(normal);
281  scalar cosAngle = max(min((cellNormal & n), 1), -1);
282  avgDiffNormal += acos(cosAngle) * mag(normal);
283  weight += mag(normal);
284  if (cosAngle < maxDiffNormal)
285  {
286  maxDiffNormal = cosAngle;
287  }
288  }
289  }
290 
291  if (weight != 0)
292  {
293  avgDiffNormal /= weight;
294  }
295  else
296  {
297  avgDiffNormal = 0;
298  }
299 
300  vector newCellNormal = normalised(interfaceNormal_[i]);
301 
302  scalar normalRes = (1 - (cellNormal & newCellNormal));
303  normalResidual[i].celli = celli;
304  normalResidual[i].normalResidual = normalRes;
305  normalResidual[i].avgAngle = avgDiffNormal;
306  }
307 }
308 
309 
310 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
311 
312 Foam::reconstruction::plicRDF::plicRDF
313 (
315  const surfaceScalarField& phi,
316  const volVectorField& U,
317  const dictionary& dict
318 )
319 :
321  (
322  typeName,
323  alpha1,
324  phi,
325  U,
326  dict
327  ),
328  mesh_(alpha1.mesh()),
329 
330  interfaceNormal_(0.2*mesh_.nCells()),
331 
332  isoFaceTol_(modelDict().getOrDefault<scalar>("isoFaceTol", 1e-8)),
333  surfCellTol_(modelDict().getOrDefault<scalar>("surfCellTol", 1e-8)),
334  tol_(modelDict().getOrDefault<scalar>("tol", 1e-6)),
335  relTol_(modelDict().getOrDefault<scalar>("relTol", 0.1)),
336  iteration_(modelDict().getOrDefault<label>("iterations", 5)),
337  interpolateNormal_(modelDict().getOrDefault("interpolateNormal", true)),
338  RDF_(mesh_),
339  sIterPLIC_(mesh_,surfCellTol_)
340 {
341  setInitNormals(false);
342 
343  centre_ = dimensionedVector("centre", dimLength, Zero);
344  normal_ = dimensionedVector("normal", dimArea, Zero);
345 
346  forAll(interfaceLabels_, i)
347  {
348  const label celli = interfaceLabels_[i];
349  if (mag(interfaceNormal_[i]) == 0)
350  {
351  continue;
352  }
353  sIterPLIC_.vofCutCell
354  (
355  celli,
356  alpha1_[celli],
357  isoFaceTol_,
358  100,
359  interfaceNormal_[i]
360  );
361 
362  if (sIterPLIC_.cellStatus() == 0)
363  {
364  normal_[celli] = sIterPLIC_.surfaceArea();
365  centre_[celli] = sIterPLIC_.surfaceCentre();
366  if (mag(normal_[celli]) == 0)
367  {
368  normal_[celli] = Zero;
369  centre_[celli] = Zero;
370  }
371  }
372  else
373  {
374  normal_[celli] = Zero;
375  centre_[celli] = Zero;
376  }
377  }
378 }
379 
380 
381 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
382 
384 {
385  addProfilingInFunction(geometricVoF);
386  zoneDistribute& exchangeFields = zoneDistribute::New(mesh_);
387  const bool uptodate = alreadyReconstructed(forceUpdate);
388 
389  if (uptodate && !forceUpdate)
390  {
391  return;
392  }
393 
394  if (mesh_.topoChanging())
395  {
396  // Introduced resizing to cope with changing meshes
397  if (interfaceCell_.size() != mesh_.nCells())
398  {
399  interfaceCell_.resize(mesh_.nCells());
400  }
401  }
402  interfaceCell_ = false;
403 
404  // Sets interfaceCell_ and interfaceNormal
405  setInitNormals(interpolateNormal_);
406 
407  centre_ = dimensionedVector("centre", dimLength, Zero);
408  normal_ = dimensionedVector("normal", dimArea, Zero);
409 
410  // nextToInterface is update on setInitNormals
411  const boolList& nextToInterface_ = RDF_.nextToInterface();
412 
413  bitSet tooCoarse(mesh_.nCells(),false);
414 
415  for (label iter=0; iter<iteration_; ++iter)
416  {
417  forAll(interfaceLabels_, i)
418  {
419  const label celli = interfaceLabels_[i];
420  if (mag(interfaceNormal_[i]) == 0 || tooCoarse.test(celli))
421  {
422  continue;
423  }
424  sIterPLIC_.vofCutCell
425  (
426  celli,
427  alpha1_[celli],
428  isoFaceTol_,
429  100,
430  interfaceNormal_[i]
431  );
432 
433  if (sIterPLIC_.cellStatus() == 0)
434  {
435 
436  normal_[celli] = sIterPLIC_.surfaceArea();
437  centre_[celli] = sIterPLIC_.surfaceCentre();
438  if (mag(normal_[celli]) == 0)
439  {
440  normal_[celli] = Zero;
441  centre_[celli] = Zero;
442  }
443  }
444  else
445  {
446  normal_[celli] = Zero;
447  centre_[celli] = Zero;
448  }
449  }
450 
451  normal_.correctBoundaryConditions();
452  centre_.correctBoundaryConditions();
453  List<normalRes> normalResidual(interfaceLabels_.size());
454 
455  surfaceVectorField::Boundary nHatb(mesh_.Sf().boundaryField());
456  nHatb *= 1/(mesh_.magSf().boundaryField());
457 
458  {
459  RDF_.constructRDF
460  (
461  nextToInterface_,
462  centre_,
463  normal_,
464  exchangeFields,
465  false
466  );
467  RDF_.updateContactAngle(alpha1_, U_, nHatb);
468  gradSurf(RDF_);
469  calcResidual(normalResidual);
470  }
471 
472  label resCounter = 0;
473  scalar avgRes = 0;
474  scalar avgNormRes = 0;
475 
476  forAll(normalResidual,i)
477  {
478 
479  const label celli = normalResidual[i].celli;
480  const scalar normalRes= normalResidual[i].normalResidual;
481  const scalar avgA = normalResidual[i].avgAngle;
482 
483  if (avgA > 0.26 && iter > 0) // 15 deg
484  {
485  tooCoarse.set(celli);
486  }
487  else
488  {
489  avgRes += normalRes;
490  scalar normRes = 0;
491  scalar discreteError = 0.01*sqr(avgA);
492  if (discreteError != 0)
493  {
494  normRes= normalRes/max(discreteError, tol_);
495  }
496  else
497  {
498  normRes= normalRes/tol_;
499  }
500  avgNormRes += normRes;
501  resCounter++;
502 
503  }
504  }
505 
506  reduce(avgRes,sumOp<scalar>());
507  reduce(avgNormRes,sumOp<scalar>());
508  reduce(resCounter,sumOp<label>());
509 
510  if (resCounter == 0) // avoid division by zero and leave loop
511  {
512  resCounter = 1;
513  avgRes = 0;
514  avgNormRes = 0;
515  }
516 
517 
518  if (iter == 0)
519  {
520  DebugInfo
521  << "initial residual absolute = "
522  << avgRes/resCounter << nl
523  << "initial residual normalized = "
524  << avgNormRes/resCounter << nl;
525  }
526 
527  if
528  (
529  (
530  (avgNormRes/resCounter < relTol_ || avgRes/resCounter < tol_)
531  && (iter >= 1 )
532  )
533  || iter + 1 == iteration_
534  )
535  {
536  DebugInfo
537  << "iterations = " << iter << nl
538  << "final residual absolute = "
539  << avgRes/resCounter << nl
540  << "final residual normalized = " << avgNormRes/resCounter
541  << endl;
542 
543  break;
544  }
545  }
546 }
547 
548 
550 {
551  // without it we seem to get a race condition
552  mesh_.C();
553 
554  cutCellPLIC cutCell(mesh_);
555 
556  forAll(normal_, celli)
557  {
558  if (mag(normal_[celli]) != 0)
559  {
560  vector n = normal_[celli]/mag(normal_[celli]);
561  scalar cutValue = (centre_[celli] - mesh_.C()[celli]) & (n);
562  cutCell.calcSubCell
563  (
564  celli,
565  cutValue,
566  n
567  );
568  alpha1_[celli] = cutCell.VolumeOfFluid();
569 
570  }
571  }
572  alpha1_.correctBoundaryConditions();
573  alpha1_.oldTime () = alpha1_;
574  alpha1_.oldTime().correctBoundaryConditions();
575 }
576 
577 
578 // ************************************************************************* //
profiling.H
Foam::reconstruction::defineTypeNameAndDebug
defineTypeNameAndDebug(isoAlpha, 0)
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::reconstructionSchemes::U_
const volVectorField & U_
Reference to the velocity field.
Definition: reconstructionSchemes.H:79
Foam::Tensor::I
static const Tensor I
Definition: Tensor.H:85
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
leastSquareGrad.H
Foam::reconstruction::plicRDF::reconstruct
virtual void reconstruct(bool forceUpdate=true)
Reconstruct interface.
Definition: plicRDF.C:383
Foam::cutCellPLIC
Class for cutting a cell, cellI, of an fvMesh, mesh_, at its intersection with an surface defined by ...
Definition: cutCellPLIC.H:71
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:65
Foam::reconstructionSchemes::normal_
volVectorField normal_
Interface area normals.
Definition: reconstructionSchemes.H:82
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::polyMesh::geometricD
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:858
Foam::reconstructionSchemes::interfaceCell_
boolList interfaceCell_
Is interface cell?
Definition: reconstructionSchemes.H:88
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::sumOp
Definition: ops.H:213
Foam::TimeState::deltaTValue
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:43
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
interpolationCellPoint.H
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::reconstructionSchemes::centre_
volVectorField centre_
Interface centres.
Definition: reconstructionSchemes.H:85
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
plicRDF.H
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::reconstructionSchemes::alpha1_
volScalarField & alpha1_
Reference to the VoF Field.
Definition: reconstructionSchemes.H:73
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:341
Foam::reconstruction::plicRDF::mapAlphaField
virtual void mapAlphaField() const
Map VoF Field in case of refinement.
Definition: plicRDF.C:549
Foam::cutCell
Service routines for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with a surface.
Definition: cutCell.H:59
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
addProfilingInFunction
#define addProfilingInFunction(name)
Definition: profilingTrigger.H:124
Foam::reconstruction::addToRunTimeSelectionTable
addToRunTimeSelectionTable(reconstructionSchemes, isoAlpha, components)
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::zoneDistribute::New
static zoneDistribute & New(const fvMesh &)
Definition: zoneDistribute.C:105
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
U
U
Definition: pEqn.H:72
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::zoneDistribute
Class for parallel communication in a narrow band. It either provides a Map with the neighbouring val...
Definition: zoneDistribute.H:64
Foam::Vector< scalar >
Foam::List< bool >
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
fvc.H
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::reconstructionSchemes
Original code supplied by Henning Scheufler, DLR (2019)
Definition: reconstructionSchemes.H:58
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::reconstructionSchemes::centre
const volVectorField & centre() const noexcept
const-Reference to interface centres
Definition: reconstructionSchemes.H:233
Foam::reconstructionSchemes::interfaceLabels_
DynamicField< label > interfaceLabels_
List of cell labels that have an interface.
Definition: reconstructionSchemes.H:91