Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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"
30#include "fvc.H"
31#include "leastSquareGrad.H"
32#include "zoneDistribute.H"
34#include "profiling.H"
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace reconstruction
41{
44}
45}
46
47
48// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49
50void Foam::reconstruction::plicRDF::interpolateNormal()
51{
52 addProfilingInFunction(geometricVoF);
53 scalar dt = mesh_.time().deltaTValue();
54
55 leastSquareGrad<scalar> lsGrad("polyDegree1",mesh_.geometricD());
56
57 zoneDistribute& exchangeFields = zoneDistribute::New(mesh_);
58 exchangeFields.setUpCommforZone(interfaceCell_, false);
59
60 Map<vector> mapCentre
61 (
63 );
64 Map<vector> mapNormal
65 (
67 );
68
69 Map<vector> mapCC
70 (
71 exchangeFields.getDatafromOtherProc(interfaceCell_, mesh_.C())
72 );
73 Map<scalar> mapAlpha
74 (
76 );
77
78 DynamicField<vector> cellCentre(100);
79 DynamicField<scalar> alphaValues(100);
80
81 DynamicList<vector> foundNormals(30);
82
83 const labelListList& stencil = exchangeFields.getStencil();
84
86 {
87 const label celli = interfaceLabels_[i];
88 vector estimatedNormal{Zero};
89 scalar weight{0};
90 foundNormals.clear();
91 for (const label gblIdx : stencil[celli])
92 {
93 vector n =
94 exchangeFields.getValue(normal_, mapNormal, gblIdx);
95 point p = mesh_.C()[celli]-U_[celli]*dt;
96 if (mag(n) != 0)
97 {
98 n /= mag(n);
100 (
101 exchangeFields.getPosition
102 (
103 centre_,
104 mapCentre,
105 gblIdx,
106 exchangeFields.getCyclicPatches
107 (
108 celli,
109 gblIdx,
110 exchangeFields.getValue(mesh_.C(), mapCC, gblIdx)
111 )
112 )
113 );
114 vector distanceToIntSeg = (tensor::I- n*n) & (p - centre);
115 estimatedNormal += n /max(mag(distanceToIntSeg), SMALL);
116 weight += 1/max(mag(distanceToIntSeg), SMALL);
117 foundNormals.append(n);
118 }
119 }
120
121 if (weight != 0 && mag(estimatedNormal) != 0)
122 {
123 estimatedNormal /= weight;
124 estimatedNormal /= mag(estimatedNormal);
125 }
126
127 bool tooCoarse = false;
128
129 if (foundNormals.size() > 1 && mag(estimatedNormal) != 0)
130 {
131 forAll(foundNormals, i)
132 {
133 // all have the length of 1
134 // to coarse if normal angle is bigger than 10 deg
135 if ((estimatedNormal & foundNormals[i]) <= 0.98)
136 {
137 tooCoarse = true;
138 }
139 }
140 }
141 else
142 {
143 tooCoarse = true;
144 }
145
146 // if a normal was found and the interface is fine enough
147 // smallDist is always smallDist
148 if (mag(estimatedNormal) != 0 && !tooCoarse)
149 {
150 interfaceNormal_[i] = estimatedNormal;
151 }
152 else
153 {
154 cellCentre.clear();
155 alphaValues.clear();
156
157 forAll(stencil[celli],i)
158 {
159 const label gblIdx = stencil[celli][i];
160 cellCentre.append
161 (
162 exchangeFields.getPosition
163 (
164 mesh_.C(), mapCC, gblIdx,
165 exchangeFields.getCyclicPatches
166 (
167 celli,
168 gblIdx,
169 exchangeFields.getValue(mesh_.C(), mapCC, gblIdx)
170 )
171 )
172 );
173 alphaValues.append
174 (
175 exchangeFields.getValue(alpha1_, mapAlpha, gblIdx)
176 );
177 }
178 cellCentre -= mesh_.C()[celli];
179 interfaceNormal_[i] = lsGrad.grad(cellCentre, alphaValues);
180 }
181
182 }
183}
184
185void Foam::reconstruction::plicRDF::gradSurf(const volScalarField& phi)
186{
187 addProfilingInFunction(geometricVoF);
188 leastSquareGrad<scalar> lsGrad("polyDegree1", mesh_.geometricD());
189 zoneDistribute& exchangeFields = zoneDistribute::New(mesh_);
190
191 exchangeFields.setUpCommforZone(interfaceCell_, false);
192
193 Map<vector> mapCC
194 (
195 exchangeFields.getDatafromOtherProc(interfaceCell_, mesh_.C())
196 );
197 Map<scalar> mapPhi
198 (
199 exchangeFields.getDatafromOtherProc(interfaceCell_, phi)
200 );
201
202 DynamicField<vector> cellCentre(100);
203 DynamicField<scalar> phiValues(100);
204
205 const labelListList& stencil = exchangeFields.getStencil();
206
207 forAll(interfaceLabels_, i)
208 {
209 const label celli = interfaceLabels_[i];
210
211 cellCentre.clear();
212 phiValues.clear();
213
214 for (const label gblIdx : stencil[celli])
215 {
216 cellCentre.append
217 (
218 exchangeFields.getPosition
219 (
220 mesh_.C(), mapCC, gblIdx,
221 exchangeFields.getCyclicPatches
222 (
223 celli,
224 gblIdx,
225 exchangeFields.getValue(mesh_.C(), mapCC, gblIdx)
226 )
227 )
228 );
229 phiValues.append
230 (
231 exchangeFields.getValue(phi, mapPhi, gblIdx)
232 );
233 }
234
235 cellCentre -= mesh_.C()[celli];
236 interfaceNormal_[i] = lsGrad.grad(cellCentre, phiValues);
237 }
238}
239
240
241void Foam::reconstruction::plicRDF::setInitNormals(bool interpolate)
242{
243 addProfilingInFunction(geometricVoF);
244 zoneDistribute& exchangeFields = zoneDistribute::New(mesh_);
245
246 interfaceLabels_.clear();
247
248 forAll(alpha1_, celli)
249 {
250 if (sIterPLIC_.isASurfaceCell(alpha1_[celli]))
251 {
252 interfaceCell_[celli] = true; // is set to false earlier
253 interfaceLabels_.append(celli);
254 }
255 }
256 interfaceNormal_.setSize(interfaceLabels_.size());
257
258 RDF_.markCellsNearSurf(interfaceCell_, 1);
259 const boolList& nextToInterface_ = RDF_.nextToInterface();
260 exchangeFields.updateStencil(nextToInterface_);
261
262 if (interpolate)
263 {
264 interpolateNormal();
265 }
266 else
267 {
268 gradSurf(alpha1_);
269 }
270}
271
272
273void Foam::reconstruction::plicRDF::calcResidual
274(
275 List<normalRes>& normalResidual
276)
277{
278 addProfilingInFunction(geometricVoF);
279 zoneDistribute& exchangeFields = zoneDistribute::New(mesh_);
280 exchangeFields.setUpCommforZone(interfaceCell_,false);
281
282 Map<vector> mapNormal
283 (
284 exchangeFields.getDatafromOtherProc(interfaceCell_, normal_)
285 );
286
287 const labelListList& stencil = exchangeFields.getStencil();
288
289 forAll(interfaceLabels_, i)
290 {
291 const label celli = interfaceLabels_[i];
292 if (mag(normal_[celli]) == 0 || mag(interfaceNormal_[i]) == 0)
293 {
294 normalResidual[i].celli = celli;
295 normalResidual[i].normalResidual = 0;
296 normalResidual[i].avgAngle = 0;
297 continue;
298 }
299
300 scalar avgDiffNormal = 0;
301 scalar maxDiffNormal = GREAT;
302 scalar weight= 0;
303 const vector cellNormal = normal_[celli]/mag(normal_[celli]);
304 forAll(stencil[celli],j)
305 {
306 const label gblIdx = stencil[celli][j];
307 vector normal =
308 exchangeFields.getValue(normal_, mapNormal, gblIdx);
309
310 if (mag(normal) != 0 && j != 0)
311 {
312 vector n = normal/mag(normal);
313 scalar cosAngle = clamp((cellNormal & n), -1, 1);
314 avgDiffNormal += acos(cosAngle) * mag(normal);
315 weight += mag(normal);
316 if (cosAngle < maxDiffNormal)
317 {
318 maxDiffNormal = cosAngle;
319 }
320 }
321 }
322
323 if (weight != 0)
324 {
325 avgDiffNormal /= weight;
326 }
327 else
328 {
329 avgDiffNormal = 0;
330 }
331
332 vector newCellNormal = normalised(interfaceNormal_[i]);
333
334 scalar normalRes = (1 - (cellNormal & newCellNormal));
335 normalResidual[i].celli = celli;
336 normalResidual[i].normalResidual = normalRes;
337 normalResidual[i].avgAngle = avgDiffNormal;
338 }
339}
340
341
342// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
343
344Foam::reconstruction::plicRDF::plicRDF
345(
347 const surfaceScalarField& phi,
348 const volVectorField& U,
349 const dictionary& dict
350)
351:
353 (
354 typeName,
355 alpha1,
356 phi,
357 U,
358 dict
359 ),
360 mesh_(alpha1.mesh()),
361
362 interfaceNormal_(0.2*mesh_.nCells()),
363
364 isoFaceTol_(modelDict().getOrDefault<scalar>("isoFaceTol", 1e-8)),
365 surfCellTol_(modelDict().getOrDefault<scalar>("surfCellTol", 1e-8)),
366 tol_(modelDict().getOrDefault<scalar>("tol", 1e-6)),
367 relTol_(modelDict().getOrDefault<scalar>("relTol", 0.1)),
368 iteration_(modelDict().getOrDefault<label>("iterations", 5)),
369 interpolateNormal_(modelDict().getOrDefault("interpolateNormal", true)),
370 RDF_(mesh_),
371 sIterPLIC_(mesh_,surfCellTol_)
372{
373 setInitNormals(false);
374
375 centre_ = Zero;
376 normal_ = Zero;
377
379 {
380 const label celli = interfaceLabels_[i];
381 if (mag(interfaceNormal_[i]) == 0)
382 {
383 continue;
384 }
385 sIterPLIC_.vofCutCell
386 (
387 celli,
388 alpha1_[celli],
389 isoFaceTol_,
390 100,
391 interfaceNormal_[i]
392 );
393
394 if (sIterPLIC_.cellStatus() == 0)
395 {
396 normal_[celli] = sIterPLIC_.surfaceArea();
397 centre_[celli] = sIterPLIC_.surfaceCentre();
398 if (mag(normal_[celli]) == 0)
399 {
400 normal_[celli] = Zero;
401 centre_[celli] = Zero;
402 }
403 }
404 else
405 {
406 normal_[celli] = Zero;
407 centre_[celli] = Zero;
409 }
410}
411
412
413// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
414
416{
417 addProfilingInFunction(geometricVoF);
418 zoneDistribute& exchangeFields = zoneDistribute::New(mesh_);
419 const bool uptodate = alreadyReconstructed(forceUpdate);
420
421 if (uptodate && !forceUpdate)
422 {
423 return;
424 }
425
426 if (mesh_.topoChanging())
427 {
428 // Introduced resizing to cope with changing meshes
429 if (interfaceCell_.size() != mesh_.nCells())
430 {
431 interfaceCell_.resize(mesh_.nCells());
432 }
433 }
434 interfaceCell_ = false;
435
436 // Sets interfaceCell_ and interfaceNormal
437 setInitNormals(interpolateNormal_);
438
439 centre_ = Zero;
440 normal_ = Zero;
441
442 // nextToInterface is update on setInitNormals
443 const boolList& nextToInterface_ = RDF_.nextToInterface();
444
445 bitSet tooCoarse(mesh_.nCells(),false);
446
447 for (label iter=0; iter<iteration_; ++iter)
448 {
449 forAll(interfaceLabels_, i)
450 {
451 const label celli = interfaceLabels_[i];
452 if (mag(interfaceNormal_[i]) == 0 || tooCoarse.test(celli))
453 {
454 continue;
455 }
456 sIterPLIC_.vofCutCell
457 (
458 celli,
459 alpha1_[celli],
460 isoFaceTol_,
461 100,
462 interfaceNormal_[i]
463 );
464
465 if (sIterPLIC_.cellStatus() == 0)
466 {
467
468 normal_[celli] = sIterPLIC_.surfaceArea();
469 centre_[celli] = sIterPLIC_.surfaceCentre();
470 if (mag(normal_[celli]) == 0)
471 {
472 normal_[celli] = Zero;
473 centre_[celli] = Zero;
474 }
475 }
476 else
477 {
478 normal_[celli] = Zero;
479 centre_[celli] = Zero;
480 }
481 }
482
483 normal_.correctBoundaryConditions();
484 centre_.correctBoundaryConditions();
485 List<normalRes> normalResidual(interfaceLabels_.size());
486
487 surfaceVectorField::Boundary nHatb(mesh_.Sf().boundaryField());
488 nHatb *= 1/(mesh_.magSf().boundaryField());
489
490 {
491 RDF_.constructRDF
492 (
493 nextToInterface_,
494 centre_,
495 normal_,
496 exchangeFields,
497 false
498 );
499 RDF_.updateContactAngle(alpha1_, U_, nHatb);
500 gradSurf(RDF_);
501 calcResidual(normalResidual);
502 }
503
504 label resCounter = 0;
505 scalar avgRes = 0;
506 scalar avgNormRes = 0;
507
508 forAll(normalResidual,i)
509 {
510
511 const label celli = normalResidual[i].celli;
512 const scalar normalRes= normalResidual[i].normalResidual;
513 const scalar avgA = normalResidual[i].avgAngle;
514
515 if (avgA > 0.26 && iter > 0) // 15 deg
516 {
517 tooCoarse.set(celli);
518 }
519 else
520 {
521 avgRes += normalRes;
522 scalar normRes = 0;
523 scalar discreteError = 0.01*sqr(avgA);
524 if (discreteError != 0)
525 {
526 normRes= normalRes/max(discreteError, tol_);
527 }
528 else
529 {
530 normRes= normalRes/tol_;
531 }
532 avgNormRes += normRes;
533 resCounter++;
534
535 }
536 }
537
538 reduce(avgRes,sumOp<scalar>());
539 reduce(avgNormRes,sumOp<scalar>());
540 reduce(resCounter,sumOp<label>());
541
542 if (resCounter == 0) // avoid division by zero and leave loop
543 {
544 resCounter = 1;
545 avgRes = 0;
546 avgNormRes = 0;
547 }
548
549
550 if (iter == 0)
551 {
553 << "initial residual absolute = "
554 << avgRes/resCounter << nl
555 << "initial residual normalized = "
556 << avgNormRes/resCounter << nl;
557 }
558
559 if
560 (
561 (
562 (avgNormRes/resCounter < relTol_ || avgRes/resCounter < tol_)
563 && (iter >= 1 )
564 )
565 || iter + 1 == iteration_
566 )
567 {
569 << "iterations = " << iter << nl
570 << "final residual absolute = "
571 << avgRes/resCounter << nl
572 << "final residual normalized = " << avgNormRes/resCounter
573 << endl;
575 break;
576 }
577 }
578}
579
580
582{
583 // without it we seem to get a race condition
584 mesh_.C();
585
586 cutCellPLIC cutCell(mesh_);
587
588 forAll(normal_, celli)
589 {
590 if (mag(normal_[celli]) != 0)
591 {
592 vector n = normal_[celli]/mag(normal_[celli]);
593 scalar cutValue = (centre_[celli] - mesh_.C()[celli]) & (n);
594 cutCell.calcSubCell
595 (
596 celli,
597 cutValue,
598 n
599 );
600 alpha1_[celli] = cutCell.VolumeOfFluid();
601
602 }
603 }
604 alpha1_.correctBoundaryConditions();
605 alpha1_.oldTime () = alpha1_;
606 alpha1_.oldTime().correctBoundaryConditions();
607}
608
609
610// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const volScalarField & alpha1
Dynamically sized Field. Similar to DynamicList, but inheriting from a Field instead of a List.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
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
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
static const Tensor I
Definition Tensor.H:81
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition bitSet.H:334
Class for cutting a cell, cellI, of an fvMesh, mesh_, at its intersection with an surface defined by ...
Definition cutCellPLIC.H:70
Service routines for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with a surface.
Definition cutCell.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
const volVectorField & C() const
Return cell centres as volVectorField.
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Estimates the gradient with a least square scheme in a cell.
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition polyMesh.C:850
Original code supplied by Henning Scheufler, DLR (2019).
boolList interfaceCell_
Is interface cell?
volVectorField normal_
Interface area normals.
const volVectorField & U_
Reference to the velocity field.
bool alreadyReconstructed(bool forceUpdate=true) const
Is the interface already up-to-date?
const volVectorField & centre() const noexcept
const-Reference to interface centres
DynamicField< label > interfaceLabels_
List of cell labels that have an interface.
reconstructionSchemes(const reconstructionSchemes &)=delete
No copy construct.
volScalarField & alpha1_
Reference to the VoF Field.
volVectorField centre_
Interface centres.
dictionary & modelDict() noexcept
Access to the model dictionary.
Reconstructs an interface (centre and normal vector) consisting of planes to match the internal fluid...
Definition plicRDF.H:72
virtual void reconstruct(bool forceUpdate=true)
Reconstruct interface.
Definition plicRDF.C:408
virtual void mapAlphaField() const
Map VoF Field in case of refinement.
Definition plicRDF.C:574
const point & surfaceCentre() const
The centre of cutting isosurface.
label cellStatus()
The cellStatus.
label vofCutCell(const label celli, const scalar alpha1, const scalar tol, const label maxIter, vector normal)
Finds matching cutValue for the given value fraction.
const vector & surfaceArea() const
The area vector of cutting isosurface.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Class for parallel communication in a narrow band. It either provides a Map with the neighbouring val...
List< label > getCyclicPatches(const label celli, const label globalIdx, const vector globalIdxCellCentre) const
Finds and returns list of all cyclic patch labels to which celli's.
void setUpCommforZone(const boolList &zone, bool updateStencil=true)
Update stencil with boolList the size has to match mesh nCells.
Type getValue(const VolumeField< Type > &phi, const Map< Type > &valuesFromOtherProc, const label gblIdx) const
Gives patchNumber and patchFaceNumber for a given Geometric volume field.
static zoneDistribute & New(const fvMesh &)
Selector.
const labelListList & getStencil() noexcept
Stencil reference.
Map< Type > getDatafromOtherProc(const boolList &zone, const VolumeField< Type > &phi)
Returns stencil and provides a Map with globalNumbering.
vector getPosition(const VolumeField< vector > &positions, const Map< vector > &valuesFromOtherProc, const label gblIdx, const List< label > cyclicPatchID=List< label >()) const
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
#define DebugInfo
Report an information message using Foam::Info.
Namespace for OpenFOAM.
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
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
Vector< scalar > vector
Definition vector.H:57
dimensionedScalar acos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define addProfilingInFunction(Name)
Define profiling trigger with specified name and description corresponding to the compiler-defined fu...
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299