Loading...
Searching...
No Matches
reconstructedDistanceFunction.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 "emptyPolyPatch.H"
30#include "processorPolyPatch.H"
31#include "syncTools.H"
32#include "unitConversion.H"
33#include "wedgePolyPatch.H"
35
36// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
37
39Foam::reconstructedDistanceFunction::coupledFacesPatch() const
40{
41 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
42
43 label nCoupled = 0;
44
45 for (const polyPatch& pp : patches)
46 {
48 {
49 nCoupled += pp.size();
50 }
51 }
52 labelList nCoupledFaces(nCoupled);
53 nCoupled = 0;
54
55 for (const polyPatch& pp : patches)
56 {
58 {
59 label facei = pp.start();
60
61 forAll(pp, i)
62 {
63 nCoupledFaces[nCoupled++] = facei++;
64 }
65 }
66 }
67
69 (
70 IndirectList<face>(mesh_.faces(), std::move(nCoupledFaces)),
71 mesh_.points()
72 );
73}
74
75
77(
78 const boolList& interfaceCells,
79 const label neiRingLevel
80)
81{
82 // performance might be improved by increasing the saving last iterations
83 // cells in a Map and loop over the map
84 if (mesh_.topoChanging())
85 {
86 // Introduced resizing to cope with changing meshes
87 if (nextToInterface_.size() != mesh_.nCells())
88 {
89 nextToInterface_.resize(mesh_.nCells());
90 }
91 coupledBoundaryPoints_ = coupledFacesPatch()().meshPoints();
92 }
93
94 const labelListList& pCells = mesh_.cellPoints();
95 const labelListList& cPoints = mesh_.pointCells();
96
97 boolList alreadyMarkedPoint(mesh_.nPoints(), false);
98 nextToInterface_ = false;
99
100 // do coupled face first
101 Map<bool> syncMap;
102
103 for (label level=0;level<=neiRingLevel;level++)
104 {
105 // parallel
106 if (level > 0)
107 {
108 forAll(coupledBoundaryPoints_, i)
109 {
110 const label pi = coupledBoundaryPoints_[i];
111 forAll(mesh_.pointCells()[pi], j)
112 {
113 const label celli = cPoints[pi][j];
114 if (cellDistLevel_[celli] == level-1)
115 {
116 syncMap.insert(pi, true);
117 break;
118 }
119 }
120 }
121
122 syncTools::syncPointMap(mesh_, syncMap, orEqOp<bool>());
123
124 // mark parallel points first
125 forAllConstIters(syncMap, iter)
126 {
127 const label pi = iter.key();
128
129 if (!alreadyMarkedPoint[pi])
130 {
131 // loop over all cells attached to the point
132 forAll(cPoints[pi], j)
133 {
134 const label pCelli = cPoints[pi][j];
135 if (cellDistLevel_[pCelli] == -1)
136 {
137 cellDistLevel_[pCelli] = level;
138 nextToInterface_[pCelli] = true;
139 }
140 }
141 }
142 alreadyMarkedPoint[pi] = true;
143 }
144 }
145
146
147 forAll(cellDistLevel_, celli)
148 {
149 if (level == 0)
150 {
151 if (interfaceCells[celli])
152 {
153 cellDistLevel_[celli] = 0;
154 nextToInterface_[celli] = true;
155 }
156 else
157 {
158 cellDistLevel_[celli] = -1;
159 }
160 }
161 else
162 {
163 if (cellDistLevel_[celli] == level-1)
164 {
165 forAll(pCells[celli], i)
166 {
167 const label pI = pCells[celli][i];
168
169 if (!alreadyMarkedPoint[pI])
170 {
171 forAll(cPoints[pI], j)
172 {
173 const label pCelli = cPoints[pI][j];
174 if (cellDistLevel_[pCelli] == -1)
175 {
176 cellDistLevel_[pCelli] = level;
177 nextToInterface_[pCelli] = true;
178 }
179 }
180 }
181 alreadyMarkedPoint[pI] = true;
182 }
183 }
184 }
186 }
187}
188
189
190// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
191
193(
194 const fvMesh& mesh
195)
196:
198 (
200 (
201 "RDF",
202 mesh.time().timeName(),
203 mesh,
204 IOobject::NO_READ,
205 IOobject::AUTO_WRITE
206 ),
207 mesh,
209 ),
210 mesh_(mesh),
211 coupledBoundaryPoints_(coupledFacesPatch()().meshPoints()),
212 cellDistLevel_
213 (
215 (
216 "cellDistLevel",
217 mesh.time().timeName(),
218 mesh,
219 IOobject::NO_READ,
220 IOobject::NO_WRITE
221 ),
222 mesh,
223 dimensionedScalar("cellDistLevel", dimless, -1)
224 ),
225 nextToInterface_(mesh.nCells(), false)
226{}
227
228
229// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
230
232(
233 const boolList& nextToInterface,
234 const volVectorField& centre,
235 const volVectorField& normal,
236 zoneDistribute& distribute,
237 bool updateStencil
238)
239{
240 volScalarField& reconDistFunc = *this;
241
242 if (nextToInterface.size() != centre.size())
243 {
245 << "size of nextToInterface: " << nextToInterface.size()
246 << "size of centre:" << centre.size()
247 << "do not match. Did the mesh change?"
248 << exit(FatalError);
249 return reconDistFunc;
250 }
251
252
253 distribute.setUpCommforZone(nextToInterface, updateStencil);
254
255 Map<vector> mapCentres =
256 distribute.getDatafromOtherProc(nextToInterface, centre);
257 Map<vector> mapNormal =
258 distribute.getDatafromOtherProc(nextToInterface, normal);
259
260 const labelListList& stencil = distribute.getStencil();
261
262
263 forAll(nextToInterface, celli)
264 {
265 if (nextToInterface[celli])
266 {
267 if (mag(normal[celli]) != 0) // interface cell
268 {
269 vector n = -normal[celli]/mag(normal[celli]);
270 scalar dist = (centre[celli] - mesh_.C()[celli]) & n;
271 reconDistFunc[celli] = dist;
272 }
273 else // nextToInterfaceCell or level == 1 cell
274 {
275 scalar averageDist = 0;
276 scalar avgWeight = 0;
277 const point p = mesh_.C()[celli];
278
279 for (const label gblIdx : stencil[celli])
280 {
281 vector n = -distribute.getValue(normal, mapNormal, gblIdx);
282 if (mag(n) != 0)
283 {
284 n /= mag(n);
285 vector c
286 (
287 distribute.getPosition
288 (
289 centre,
290 mapCentres,
291 gblIdx,
292 distribute.getCyclicPatches
293 (
294 celli,
295 gblIdx,
296 distribute.getValue
297 (
298 centre,
299 mapCentres,
300 gblIdx
301 )
302 )
303 )
304 );
305 vector distanceToIntSeg(c - p);
306 scalar distToSurf = distanceToIntSeg & n;
307 scalar weight = 0;
308
309 if (mag(distanceToIntSeg) != 0)
310 {
311 distanceToIntSeg /= mag(distanceToIntSeg);
312 weight = sqr(mag(distanceToIntSeg & n));
313 }
314 else // exactly on the center
315 {
316 weight = 1;
317 }
318 averageDist += distToSurf * weight;
319 avgWeight += weight;
320 }
321 }
322
323 if (avgWeight != 0)
324 {
325 reconDistFunc[celli] = averageDist / avgWeight;
326 }
327 }
328 }
329 else
330 {
331 reconDistFunc[celli] = 0;
332 }
333 }
334
335 forAll(reconDistFunc.boundaryField(), patchI)
336 {
337 fvPatchScalarField& pRDF = reconDistFunc.boundaryFieldRef()[patchI];
339 {
340 const polyPatch& pp = pRDF.patch().patch();
341 forAll(pRDF, i)
342 {
343 const label pCellI = pp.faceCells()[i];
344
345 if (nextToInterface_[pCellI])
346 {
347 scalar averageDist = 0;
348 scalar avgWeight = 0;
349 const point p = mesh_.C().boundaryField()[patchI][i];
350
351 forAll(stencil[pCellI], j)
352 {
353 const label gblIdx = stencil[pCellI][j];
354 vector n = -distribute.getValue
355 (
356 normal,
357 mapNormal,
358 gblIdx
359 );
360 if (mag(n) != 0)
361 {
362 n /= mag(n);
363 vector c
364 (
365 distribute.getPosition
366 (
367 centre,
368 mapCentres,
369 gblIdx,
370 distribute.getCyclicPatches
371 (
372 pCellI,
373 gblIdx,
374 distribute.getValue
375 (
376 centre,
377 mapCentres,
378 gblIdx
379 )
380 )
381 )
382 );
383 vector distanceToIntSeg(c - p);
384 scalar distToSurf = distanceToIntSeg & n;
385 scalar weight = 0;
386
387 if (mag(distanceToIntSeg) != 0)
388 {
389 distanceToIntSeg /= mag(distanceToIntSeg);
390 weight = sqr(mag(distanceToIntSeg & n));
391 }
392 else // exactly on the center
393 {
394 weight = 1;
395 }
396 averageDist += distToSurf * weight;
397 avgWeight += weight;
398 }
399 }
400
401 if (avgWeight != 0)
402 {
403 pRDF[i] = averageDist / avgWeight;
404 }
405 else
406 {
407 pRDF[i] = 0;
408 }
409 }
410 else
411 {
412 pRDF[i] = 0;
413 }
414 }
415 }
416 }
418 reconDistFunc.correctBoundaryConditions();
419
420 return reconDistFunc;
421}
422
423
425(
426 const volScalarField& alpha,
427 const volVectorField& U,
429)
430{
431 const fvMesh& mesh = alpha.mesh();
432 const volScalarField::Boundary& abf = alpha.boundaryField();
434
435 const fvBoundaryMesh& boundary = mesh.boundary();
436
437 forAll(boundary, patchi)
438 {
440 {
443 (
445 (
446 abf[patchi]
447 )
448 );
449
450 fvsPatchVectorField& nHatp = nHatb[patchi];
451 const scalarField theta
452 (
453 degToRad()*acap.theta(U.boundaryField()[patchi], nHatp)
454 );
455
456 RDFbf[patchi] =
457 1/acap.patch().deltaCoeffs()*cos(theta)
458 + RDFbf[patchi].patchInternalField();
459 }
460 }
461}
462
463
464// ************************************************************************* //
constexpr scalar pi(M_PI)
label n
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const Mesh & mesh() const noexcept
Return const reference to mesh.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
void correctBoundaryConditions()
Correct boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
A List with indirect addressing.
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Abstract base class for two-phase alphaContactAngle boundary conditions.
virtual tmp< scalarField > theta(const fvPatchVectorField &Up, const fvsPatchVectorField &nHat) const =0
Return the contact angle.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvPatch & patch() const noexcept
Return the patch.
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition fvPatch.H:202
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
void markCellsNearSurf(const boolList &interfaceCells, const label neiRingLevel)
const boolList & nextToInterface() const noexcept
reconstructedDistanceFunction(const fvMesh &mesh)
Construct from fvMesh.
const volScalarField & constructRDF(const boolList &nextToInterface, const volVectorField &centre, const volVectorField &normal, zoneDistribute &distribute, bool updateStencil=true)
void updateContactAngle(const volScalarField &alpha, const volVectorField &U, surfaceVectorField::Boundary &nHatb)
static void syncPointMap(const polyMesh &mesh, Map< T > &pointValues, const CombineOp &cop, const TransformOp &top)
Synchronize values on selected points.
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.
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
U
Definition pEqn.H:72
volScalarField & p
const polyBoundaryMesh & patches
faceListList boundary
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
word timeName
Definition getTimeIndex.H:3
const dimensionedScalar c
Speed of light in a vacuum.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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
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.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
fvsPatchField< vector > fvsPatchVectorField
dimensionedScalar cos(const dimensionedScalar &ds)
fvPatchField< scalar > fvPatchScalarField
volScalarField & alpha
Us boundaryFieldRef().evaluateCoupled< coupledFaPatch >()
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Unit conversion functions.