conformalVoronoiMeshI.H
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) 2012-2015 OpenFOAM Foundation
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 "indexedVertexOps.H"
29 #include "indexedCellOps.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 inline Foam::scalar Foam::conformalVoronoiMesh::defaultCellSize() const
34 {
36 }
37 
38 
39 inline Foam::scalar Foam::conformalVoronoiMesh::targetCellSize
40 (
41  const Foam::point& pt
42 ) const
43 {
44  return cellShapeControls().cellSize(pt);
45 }
46 
47 
48 inline Foam::scalar Foam::conformalVoronoiMesh::averageAnyCellSize
49 (
50  const Vertex_handle& vA,
51  const Vertex_handle& vB
52 ) const
53 {
54  if
55  (
56  (!vA->internalOrBoundaryPoint() || vA->referred())
57  && (!vB->internalOrBoundaryPoint() || vB->referred())
58  )
59  {
60  // There are no internalOrBoundaryPoints available, determine
61  // size from scratch
62 
63  // Geometric mean
64  return sqrt
65  (
66  targetCellSize(topoint(vA->point()))
67  *targetCellSize(topoint(vB->point()))
68  );
69  }
70  else if (!vB->internalOrBoundaryPoint() || vB->referred())
71  {
72  return vA->targetCellSize();
73  }
74  else if (!vA->internalOrBoundaryPoint() || vA->referred())
75  {
76  return vB->targetCellSize();
77  }
78 
80 }
81 
82 
83 inline Foam::scalar Foam::conformalVoronoiMesh::averageAnyCellSize
84 (
85  const Delaunay::Finite_facets_iterator& fit
86 ) const
87 {
88  // Arithmetic mean
89 
90  scalar sizeSum = 0;
91 
92  label nProducts = 0;
93 
94  const Cell_handle c(fit->first);
95  const label oppositeVertex = fit->second;
96 
97  for (label i = 0; i < 3; i++)
98  {
99  Vertex_handle v = c->vertex(vertex_triple_index(oppositeVertex, i));
100 
101  if (v->internalOrBoundaryPoint() && !v->referred())
102  {
103 
104  sizeSum += v->targetCellSize();
105 
106  nProducts++;
107  }
108  }
109 
110  if (nProducts < 1)
111  {
112  // There are no internalOrBoundaryPoints available, determine
113  // size from scratch
114 
115  for (label i = 0; i < 3; i++)
116  {
117  Vertex_handle v = c->vertex(vertex_triple_index(oppositeVertex, i));
118 
119  sizeSum += targetCellSize(topoint(v->point()));
120  }
121 
122  nProducts = 3;
123  }
124 
125  if (sizeSum < 0)
126  {
128  << "sizeSum = " << sizeSum
129  << endl;
130 
131  return 0;
132  }
133 
134  return pow(sizeSum, (1.0/nProducts));
135 }
136 
137 
139 (
140  const Foam::point& pt
141 ) const
142 {
143  return targetCellSize(pt)*foamyHexMeshControls().pointPairDistanceCoeff();
144 }
145 
146 
148 (
149  const Foam::point& pt
150 ) const
151 {
152  return
153  pointPairDistance(pt)
154  *foamyHexMeshControls().mixedFeaturePointPPDistanceCoeff();
155 }
156 
157 
159 (
160  const Foam::point& pt
161 ) const
162 {
163  return
164  sqr
165  (
166  targetCellSize(pt)
167  *foamyHexMeshControls().featurePointExclusionDistanceCoeff()
168  );
169 }
170 
171 
173 (
174  const Foam::point& pt
175 ) const
176 {
177  return
178  sqr
179  (
180  targetCellSize(pt)
181  *foamyHexMeshControls().featureEdgeExclusionDistanceCoeff()
182  );
183 }
184 
185 
187 (
188  const Foam::point& pt
189 ) const
190 {
191  return
192  sqr
193  (
194  targetCellSize(pt)
195  *foamyHexMeshControls().surfacePtExclusionDistanceCoeff()
196  );
197 }
198 
199 
201 (
202  const Foam::point& pt
203 ) const
204 {
205  return
206  sqr
207  (
208  targetCellSize(pt)
209  *foamyHexMeshControls().surfaceSearchDistanceCoeff()
210  );
211 }
212 
213 
215 (
216  const Foam::point& pt
217 ) const
218 {
219  return
220  targetCellSize(pt)
221  *foamyHexMeshControls().maxSurfaceProtrusionCoeff();
222 }
223 
224 
225 inline void Foam::conformalVoronoiMesh::createPointPair
226 (
227  const scalar ppDist,
228  const Foam::point& surfPt,
229  const vector& n,
230  const bool ptPair,
231  DynamicList<Vb>& pts
232 ) const
233 {
234  vector ppDistn = ppDist*n;
235 
236 // const Foam::point internalPt = surfPt - ppDistn;
237 // const Foam::point externalPt = surfPt + ppDistn;
238 
239 // bool internalInside = geometryToConformTo_.inside(internalPt);
240 // bool externalOutside = geometryToConformTo_.outside(externalPt);
241 
242 // if (internalInside && externalOutside)
243  {
244  pts.append
245  (
246  Vb
247  (
248  surfPt - ppDistn,
249  vertexCount() + pts.size(),
250  Vb::vtInternalSurface,
252  )
253  );
254 
255  pts.append
256  (
257  Vb
258  (
259  surfPt + ppDistn,
260  vertexCount() + pts.size(),
261  Vb::vtExternalSurface,
263  )
264  );
265 
266  if (ptPair)
267  {
268  ptPairs_.addPointPair
269  (
270  pts[pts.size() - 2].index(),
271  pts[pts.size() - 1].index() // external 0 -> slave
272  );
273  }
274  }
275 // else
276 // {
277 // Info<< "Warning: point pair not inside/outside" << nl
278 // << " surfPt = " << surfPt << nl
279 // << " internal = " << internalPt << " " << internalInside << nl
280 // << " external = " << externalPt << " " << externalOutside
281 // << endl;
282 // }
283 }
284 
285 
286 inline Foam::point Foam::conformalVoronoiMesh::perturbPoint
287 (
288  const Foam::point& pt
289 ) const
290 {
291  Foam::point perturbedPt(pt);
292 
293 // vector delta(xR/ni, yR/nj, zR/nk);
294 // scalar pert = randomPerturbationCoeff*cmptMin(delta);
295 
296  scalar pert = 1e-12*defaultCellSize();
297 
298  perturbedPt.x() += pert*(rndGen_.sample01<scalar>() - 0.5);
299  perturbedPt.y() += pert*(rndGen_.sample01<scalar>() - 0.5);
300  perturbedPt.z() += pert*(rndGen_.sample01<scalar>() - 0.5);
301 
302  return perturbedPt;
303 }
304 
305 
306 inline void Foam::conformalVoronoiMesh::createBafflePointPair
307 (
308  const scalar ppDist,
309  const Foam::point& surfPt,
310  const vector& n,
311  const bool ptPair,
312  DynamicList<Vb>& pts
313 ) const
314 {
315  vector ppDistn = ppDist*n;
316 
317  pts.append
318  (
319  Vb
320  (
321  surfPt - ppDistn,
322  vertexCount() + pts.size(),
323  Vb::vtInternalSurfaceBaffle,
325  )
326  );
327 
328  pts.append
329  (
330  Vb
331  (
332  surfPt + ppDistn,
333  vertexCount() + pts.size(),
334  Vb::vtExternalSurfaceBaffle,
336  )
337  );
338 
339  if (ptPair)
340  {
341  ptPairs_.addPointPair
342  (
343  pts[pts.size() - 2].index(), // external 0 -> slave
344  pts[pts.size() - 1].index()
345  );
346  }
347 }
348 
349 
350 inline bool Foam::conformalVoronoiMesh::internalPointIsInside
351 (
352  const Foam::point& pt
353 ) const
354 {
355  if
356  (
357  !geometryToConformTo_.globalBounds().contains(pt)
358  || !geometryToConformTo_.inside(pt)
359  )
360  {
361  return false;
362  }
363 
364  return true;
365 }
366 
367 
368 inline bool Foam::conformalVoronoiMesh::isBoundaryDualFace
369 (
370  const Delaunay::Finite_edges_iterator& eit
371 ) const
372 {
373  Cell_handle c = eit->first;
374  Vertex_handle vA = c->vertex(eit->second);
375  Vertex_handle vB = c->vertex(eit->third);
376 
377 // if (vA->internalBoundaryPoint() && vB->externalBoundaryPoint())
378 // {
379 // if (vA->index() == vB->index() - 1)
380 // {
381 // return true;
382 // }
383 // }
384 // else if (vA->externalBoundaryPoint() && vB->internalBoundaryPoint())
385 // {
386 // if (vA->index() == vB->index() + 1)
387 // {
388 // return true;
389 // }
390 // }
391 //
392 // return false;
393 
394  // A dual face on the boundary will result from one Dv inside and
395  // one outside
396  return
397  (
398  (
399  (vA->internalOrBoundaryPoint() && !vA->referred())
400  || (vB->internalOrBoundaryPoint() && !vB->referred())
401  )
402  && (
403  !vA->internalOrBoundaryPoint()
404  || !vB->internalOrBoundaryPoint()
405  )
406  );
407 }
408 
409 
410 inline Foam::List<bool> Foam::conformalVoronoiMesh::dualFaceBoundaryPoints
411 (
412  const Delaunay::Finite_edges_iterator& eit
413 ) const
414 {
415  Cell_circulator ccStart = incident_cells(*eit);
416  Cell_circulator cc1 = ccStart;
417  Cell_circulator cc2 = cc1;
418 
419  // Advance the second circulator so that it always stays on the next
420  // cell around the edge;
421  cc2++;
422 
423  DynamicList<bool> tmpFaceBoundaryPoints;
424 
425  do
426  {
427  label cc1I = cc1->cellIndex();
428 
429  label cc2I = cc2->cellIndex();
430 
431  if (cc1I != cc2I)
432  {
433  if (cc1->boundaryDualVertex())
434  {
435  tmpFaceBoundaryPoints.append(true);
436  }
437  else
438  {
439  tmpFaceBoundaryPoints.append(false);
440  }
441  }
442 
443  cc1++;
444 
445  cc2++;
446 
447  } while (cc1 != ccStart);
448 
449  return tmpFaceBoundaryPoints;
450 }
451 
452 
453 inline Foam::List<Foam::label> Foam::conformalVoronoiMesh::processorsAttached
454 (
455  const Delaunay::Finite_facets_iterator& fit
456 ) const
457 {
458  DynamicList<label> procsAttached(8);
459 
460  const Cell_handle c1(fit->first);
461  const label oppositeVertex = fit->second;
462  const Cell_handle c2(c1->neighbor(oppositeVertex));
463 
464  FixedList<label, 4> c1Procs(CGAL::indexedCellOps::processorsAttached(c1));
465 
466  FixedList<label, 4> c2Procs(CGAL::indexedCellOps::processorsAttached(c2));
467 
468  forAll(c1Procs, aPI)
469  {
470  if (!procsAttached.found(c1Procs[aPI]))
471  {
472  procsAttached.append(c1Procs[aPI]);
473  }
474 
475  if (!procsAttached.found(c2Procs[aPI]))
476  {
477  procsAttached.append(c2Procs[aPI]);
478  }
479  }
480 
481  return List<label>(procsAttached);
482 }
483 
484 
485 inline bool Foam::conformalVoronoiMesh::isParallelDualEdge
486 (
487  const Delaunay::Finite_facets_iterator& fit
488 ) const
489 {
490  const Cell_handle c1(fit->first);
491  const label oppositeVertex = fit->second;
492 
493  return
494  (
495  c1->vertex(vertex_triple_index(oppositeVertex, 0))->referred()
496  || c1->vertex(vertex_triple_index(oppositeVertex, 1))->referred()
497  || c1->vertex(vertex_triple_index(oppositeVertex, 2))->referred()
498  );
499 }
500 
501 
502 inline bool Foam::conformalVoronoiMesh::isProcBoundaryEdge
503 (
504  const Delaunay::Finite_edges_iterator& eit
505 ) const
506 {
507  bool isProcBoundaryEdge = false;
508 
509  Cell_handle c = eit->first;
510  Vertex_handle vA = c->vertex(eit->second);
511  Vertex_handle vB = c->vertex(eit->third);
512 
513  if
514  (
515  (
516  (vA->referred() && !vB->referred())
517  || (vB->referred() && !vA->referred())
518  )
519  && vA->internalOrBoundaryPoint()
520  && vB->internalOrBoundaryPoint()
521  )
522  {
523  isProcBoundaryEdge = true;
524  }
525 
526  return isProcBoundaryEdge;
527 }
528 
529 
530 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
531 
533 {
534  return runTime_;
535 }
536 
537 
539 {
540  return rndGen_;
541 }
542 
543 
544 inline const Foam::searchableSurfaces&
546 {
547  return allGeometry_;
548 }
549 
550 
551 inline const Foam::conformationSurfaces&
553 {
554  return geometryToConformTo_;
555 }
556 
557 
560 {
561  if (!Pstream::parRun())
562  {
564  << "The backgroundMeshDecomposition cannot be asked for in serial."
565  << exit(FatalError) << endl;
566  }
567 
568  return *decomposition_;
569 }
570 
571 
572 inline const Foam::cellShapeControl&
574 {
575  return cellShapeControl_;
576 }
577 
578 
579 inline const Foam::cvControls&
581 {
582  return foamyHexMeshControls_;
583 }
584 
585 
586 // ************************************************************************* //
Foam::conformalVoronoiMesh::cellShapeControls
const cellShapeControl & cellShapeControls() const
Return the cellShapeControl object.
Definition: conformalVoronoiMeshI.H:573
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::cvControls
Controls for the conformalVoronoiMesh mesh generator.
Definition: cvControls.H:55
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
Foam::conformalVoronoiMesh::surfaceSearchDistanceSqr
scalar surfaceSearchDistanceSqr(const Foam::point &pt) const
Return the square of the local surface search distance.
Definition: conformalVoronoiMeshI.H:201
Foam::conformalVoronoiMesh::pointPairDistance
scalar pointPairDistance(const Foam::point &pt) const
Return the local point pair separation at the given location.
Definition: conformalVoronoiMeshI.H:139
Foam::UPstream::parRun
static bool & parRun()
Test if this a parallel run, or allow modify access.
Definition: UPstream.H:434
CGAL::indexedVertex
An indexed form of CGAL::Triangulation_vertex_base_3<K> used to keep track of the Delaunay vertices i...
Definition: indexedVertex.H:54
Foam::conformalVoronoiMesh::geometryToConformTo
const conformationSurfaces & geometryToConformTo() const
Return the conformationSurfaces object.
Definition: conformalVoronoiMeshI.H:552
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::conformalVoronoiMesh::mixedFeaturePointDistance
scalar mixedFeaturePointDistance(const Foam::point &pt) const
Return the local mixed feature point placement distance.
Definition: conformalVoronoiMeshI.H:148
Foam::conformalVoronoiMesh::rndGen
Random & rndGen() const
Return the random number generator.
Definition: conformalVoronoiMeshI.H:538
Foam::conformationSurfaces
Definition: conformationSurfaces.H:55
Foam::topoint
pointFromPoint topoint(const Point &P)
Definition: pointConversion.H:72
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::constant::physicoChemical::c1
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Foam::conformalVoronoiMesh::decomposition
const backgroundMeshDecomposition & decomposition() const
Return the backgroundMeshDecomposition.
Definition: conformalVoronoiMeshI.H:559
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:474
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::FatalError
error FatalError
Foam::cellShapeControl
Definition: cellShapeControl.H:64
Foam::constant::physicoChemical::c2
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::conformalVoronoiMesh::foamyHexMeshControls
const cvControls & foamyHexMeshControls() const
Return the foamyHexMeshControls object.
Definition: conformalVoronoiMeshI.H:580
Foam::conformalVoronoiMesh::time
const Time & time() const
Return the Time object.
Definition: conformalVoronoiMeshI.H:532
Foam::conformalVoronoiMesh::surfacePtExclusionDistanceSqr
scalar surfacePtExclusionDistanceSqr(const Foam::point &pt) const
Return the square of the local surface point exclusion distance.
Definition: conformalVoronoiMeshI.H:187
Foam::conformalVoronoiMesh::featureEdgeExclusionDistanceSqr
scalar featureEdgeExclusionDistanceSqr(const Foam::point &pt) const
Return the square of the local feature edge exclusion distance.
Definition: conformalVoronoiMeshI.H:173
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:464
Foam::Vector
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:62
Foam::List< bool >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::searchableSurfaces
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
Definition: searchableSurfaces.H:92
Foam::backgroundMeshDecomposition
Store a background polyMesh to use for the decomposition of space and queries for parallel conformalV...
Definition: backgroundMeshDecomposition.H:94
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::conformalVoronoiMesh::maxSurfaceProtrusion
scalar maxSurfaceProtrusion(const Foam::point &pt) const
Return the local maximum surface protrusion distance.
Definition: conformalVoronoiMeshI.H:215
Foam::conformalVoronoiMesh::featurePointExclusionDistanceSqr
scalar featurePointExclusionDistanceSqr(const Foam::point &pt) const
Return the square of the local feature point exclusion distance.
Definition: conformalVoronoiMeshI.H:159
CGAL::indexedVertexOps::averageCellSize
Foam::scalar averageCellSize(const VertexType &vA, const VertexType &vB)
Return the target cell size from that stored on a pair of Delaunay vertices,.
CGAL::indexedCellOps::processorsAttached
Foam::FixedList< Foam::label, 4 > processorsAttached(const CellType &c)
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::cvControls::defaultCellSize
scalar defaultCellSize() const
Return the defaultCellSize.
Definition: cvControlsI.H:139
indexedCellOps.H
Foam::conformalVoronoiMesh::allGeometry
const searchableSurfaces & allGeometry() const
Return the allGeometry object.
Definition: conformalVoronoiMeshI.H:545
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
indexedVertexOps.H