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 -------------------------------------------------------------------------------
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 "emptyPolyPatch.H"
30 #include "processorPolyPatch.H"
31 #include "syncTools.H"
32 #include "unitConversion.H"
33 #include "wedgePolyPatch.H"
35 
36 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
37 
39 Foam::reconstructedDistanceFunction::coupledFacesPatch() const
40 {
41  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
42 
43  label nCoupled = 0;
44 
45  for (const polyPatch& pp : patches)
46  {
47  if (isA<coupledPolyPatch>(pp))
48  {
49  nCoupled += pp.size();
50  }
51  }
52  labelList nCoupledFaces(nCoupled);
53  nCoupled = 0;
54 
55  for (const polyPatch& pp : patches)
56  {
57  if (isA<coupledPolyPatch>(pp))
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  }
185  }
186  }
187 }
188 
189 
190 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
191 
193 (
194  const fvMesh& mesh
195 )
196 :
198  (
199  IOobject
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  (
214  IOobject
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 = distribute.getValue(centre,mapCentres,gblIdx);
286  vector distanceToIntSeg = (c - p);
287  scalar distToSurf = distanceToIntSeg & (n);
288  scalar weight = 0;
289 
290  if (mag(distanceToIntSeg) != 0)
291  {
292  distanceToIntSeg /= mag(distanceToIntSeg);
293  weight = sqr(mag(distanceToIntSeg & n));
294  }
295  else // exactly on the center
296  {
297  weight = 1;
298  }
299  averageDist += distToSurf * weight;
300  avgWeight += weight;
301  }
302  }
303 
304  if (avgWeight != 0)
305  {
306  reconDistFunc[celli] = averageDist / avgWeight;
307  }
308  }
309  }
310  else
311  {
312  reconDistFunc[celli] = 0;
313  }
314  }
315 
316  forAll(reconDistFunc.boundaryField(), patchI)
317  {
318  fvPatchScalarField& pRDF = reconDistFunc.boundaryFieldRef()[patchI];
319  if (isA<calculatedFvPatchScalarField>(pRDF))
320  {
321  const polyPatch& pp = pRDF.patch().patch();
322  forAll(pRDF, i)
323  {
324  const label pCellI = pp.faceCells()[i];
325 
326  if (nextToInterface_[pCellI])
327  {
328  scalar averageDist = 0;
329  scalar avgWeight = 0;
330  const point p = mesh_.C().boundaryField()[patchI][i];
331 
332  forAll(stencil[pCellI], j)
333  {
334  const label gblIdx = stencil[pCellI][j];
335  vector n = -distribute.getValue(normal, mapNormal, gblIdx);
336  if (mag(n) != 0)
337  {
338  n /= mag(n);
339  vector c =
340  distribute.getValue(centre, mapCentres, gblIdx);
341  vector distanceToIntSeg = (c - p);
342  scalar distToSurf = distanceToIntSeg & (n);
343  scalar weight = 0;
344 
345  if (mag(distanceToIntSeg) != 0)
346  {
347  distanceToIntSeg /= mag(distanceToIntSeg);
348  weight = sqr(mag(distanceToIntSeg & n));
349  }
350  else // exactly on the center
351  {
352  weight = 1;
353  }
354  averageDist += distToSurf * weight;
355  avgWeight += weight;
356  }
357  }
358 
359  if (avgWeight != 0)
360  {
361  pRDF[i] = averageDist / avgWeight;
362  }
363  else
364  {
365  pRDF[i] = 0;
366  }
367  }
368  else
369  {
370  pRDF[i] = 0;
371  }
372  }
373  }
374  }
375 
376  reconDistFunc.correctBoundaryConditions();
377 
378  return reconDistFunc;
379 }
380 
381 
383 (
384  const volScalarField& alpha,
385  const volVectorField& U,
387 )
388 {
389  const fvMesh& mesh = alpha.mesh();
390  const volScalarField::Boundary& abf = alpha.boundaryField();
392 
394 
395  forAll(boundary, patchi)
396  {
397  if (isA<alphaContactAngleTwoPhaseFvPatchScalarField>(abf[patchi]))
398  {
401  (
402  refCast<const alphaContactAngleTwoPhaseFvPatchScalarField>
403  (
404  abf[patchi]
405  )
406  );
407 
408  fvsPatchVectorField& nHatp = nHatb[patchi];
409  const scalarField theta
410  (
411  degToRad()*acap.theta(U.boundaryField()[patchi], nHatp)
412  );
413 
414  RDFbf[patchi] =
415  1/acap.patch().deltaCoeffs()*cos(theta)
416  + RDFbf[patchi].patchInternalField();
417  }
418  }
419 }
420 
421 
422 // ************************************************************************* //
const labelListList & getStencil() noexcept
Stencil reference.
faceListList boundary
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
Map< Type > getDatafromOtherProc(const boolList &zone, const VolumeField< Type > &phi)
Returns stencil and provides a Map with globalNumbering.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
reconstructedDistanceFunction(const fvMesh &mesh)
Construct from fvMesh.
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
void updateContactAngle(const volScalarField &alpha, const volVectorField &U, surfaceVectorField::Boundary &nHatb)
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
fvPatchField< scalar > fvPatchScalarField
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
void markCellsNearSurf(const boolList &interfaceCells, const label neiRingLevel)
constexpr scalar pi(M_PI)
Vector< scalar > vector
Definition: vector.H:57
const volScalarField & constructRDF(const boolList &nextToInterface, const volVectorField &centre, const volVectorField &normal, zoneDistribute &distribute, bool updateStencil=true)
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
Abstract base class for two-phase alphaContactAngle boundary conditions.
Us boundaryFieldRef().evaluateCoupled< coupledFaPatch >()
void setUpCommforZone(const boolList &zone, bool updateStencil=true)
Update stencil with boolList the size has to match mesh nCells.
virtual tmp< scalarField > theta(const fvPatchVectorField &Up, const fvsPatchVectorField &nHat) const =0
Return the contact angle.
U
Definition: pEqn.H:72
vector point
Point is a vector.
Definition: point.H:37
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Foam::fvBoundaryMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionedScalar c
Speed of light in a vacuum.
const polyBoundaryMesh & patches
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
label n
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Class for parallel communication in a narrow band. It either provides a Map with the neighbouring val...
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
List< bool > boolList
A List of bools.
Definition: List.H:60
A List with indirect addressing.
Definition: IndirectList.H:60
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
static void syncPointMap(const polyMesh &mesh, Map< T > &pointValues, const CombineOp &cop, const TransformOp &top)
Synchronize values on selected points.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
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 constexpr const zero Zero
Global zero (0)
Definition: zero.H:127