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