conformalVoronoiMeshConformToSurface.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) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2020-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "conformalVoronoiMesh.H"
31 #include "vectorTools.H"
32 #include "indexedCellChecks.H"
33 #include "IOmanip.H"
34 #include "OBJstream.H"
35 
36 using namespace Foam::vectorTools;
37 
38 const Foam::scalar Foam::conformalVoronoiMesh::searchConeAngle
39  = Foam::cos(degToRad(30.0));
40 
41 const Foam::scalar Foam::conformalVoronoiMesh::searchAngleOppositeSurface
42  = Foam::cos(degToRad(150.0));
43 
44 
45 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
46 
47 void Foam::conformalVoronoiMesh::conformToSurface()
48 {
49  this->resetCellCount();
50  // Index the cells
51  for
52  (
53  Delaunay::Finite_cells_iterator cit = finite_cells_begin();
54  cit != finite_cells_end();
55  ++cit
56  )
57  {
58  cit->cellIndex() = Cb::ctUnassigned;
59  }
60 
61  if (!reconformToSurface())
62  {
63  // Reinsert stored surface conformation
64  reinsertSurfaceConformation();
65 
66  if (Pstream::parRun())
67  {
68  sync(decomposition().procBounds());
69  }
70  }
71  else
72  {
73  ptPairs_.clear();
74 
75  // Rebuild, insert and store new surface conformation
76  buildSurfaceConformation();
77 
78  if (distributeBackground(*this))
79  {
80  if (Pstream::parRun())
81  {
82  sync(decomposition().procBounds());
83  }
84  }
85 
86  // Do not store the surface conformation until after it has been
87  // (potentially) redistributed.
88  storeSurfaceConformation();
89  }
90 
91  // reportSurfaceConformationQuality();
92 }
93 
94 
95 bool Foam::conformalVoronoiMesh::reconformToSurface() const
96 {
97  if
98  (
99  runTime_.timeIndex()
100  % foamyHexMeshControls().surfaceConformationRebuildFrequency() == 0
101  )
102  {
103  return true;
104  }
105 
106  return false;
107 }
108 
109 
110 // TODO: Investigate topological tests
111 Foam::label Foam::conformalVoronoiMesh::findVerticesNearBoundaries()
112 {
113  label countNearBoundaryVertices = 0;
114 
115  for
116  (
117  Delaunay::Finite_facets_iterator fit = finite_facets_begin();
118  fit != finite_facets_end();
119  ++fit
120  )
121  {
122  Cell_handle c1 = fit->first;
123  Cell_handle c2 = fit->first->neighbor(fit->second);
124 
125  if (is_infinite(c1) || is_infinite(c2))
126  {
127  continue;
128  }
129 
130  pointFromPoint dE0 = c1->dual();
131  pointFromPoint dE1 = c2->dual();
132 
133  if (!geometryToConformTo_.findSurfaceAnyIntersection(dE0, dE1))
134  {
135  continue;
136  }
137 
138  for (label celli = 0; celli < 4; ++celli)
139  {
140  Vertex_handle v = c1->vertex(celli);
141 
142  if
143  (
144  !is_infinite(v)
145  && v->internalPoint()
146  && fit->second != celli
147  )
148  {
149  v->setNearBoundary();
150  }
151  }
152 
153  for (label celli = 0; celli < 4; ++celli)
154  {
155  Vertex_handle v = c2->vertex(celli);
156 
157  if
158  (
159  !is_infinite(v)
160  && v->internalPoint()
161  && fit->second != celli
162  )
163  {
164  v->setNearBoundary();
165  }
166  }
167  }
168 
169  for
170  (
171  Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
172  vit != finite_vertices_end();
173  ++vit
174  )
175  {
176  if (vit->nearBoundary())
177  {
178  countNearBoundaryVertices++;
179  }
180  }
181 
182  // Geometric test.
183 // for
184 // (
185 // Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
186 // vit != finite_vertices_end();
187 // ++vit
188 // )
189 // {
190 // if (vit->internalPoint() && !vit->nearBoundary())
191 // {
192 // pointFromPoint pt = topoint(vit->point());
193 //
194 // const scalar range = sqr
195 // (
196 // foamyHexMeshControls().nearBoundaryDistanceCoeff()
197 // *targetCellSize(pt)
198 // );
199 //
200 // pointIndexHit pHit;
201 // label hitSurface;
202 //
203 // geometryToConformTo_.findSurfaceNearest
204 // (
205 // pt,
206 // range,
207 // pHit,
208 // hitSurface
209 // );
210 //
211 // if (pHit.hit())
212 // {
213 // vit->setNearBoundary();
214 // countNearBoundaryVertices++;
215 // }
216 // }
217 // }
218 
219  return countNearBoundaryVertices;
220 }
221 
222 
223 void Foam::conformalVoronoiMesh::buildSurfaceConformation()
224 {
225  timeCheck("Start buildSurfaceConformation");
226 
227  Info<< nl
228  << "Rebuilding surface conformation for more iterations"
229  << endl;
230 
231  existingEdgeLocations_.clearStorage();
232  existingSurfacePtLocations_.clearStorage();
233 
234  buildEdgeLocationTree(existingEdgeLocations_);
235  buildSurfacePtLocationTree(existingSurfacePtLocations_);
236 
237  label initialTotalHits = 0;
238 
239  // Surface protrusion conformation is done in two steps.
240  // 1. the dual edges (of all internal vertices) can stretch to
241  // 'infinity' so any intersection would be badly behaved. So
242  // just find the nearest point on the geometry and insert point
243  // pairs.
244  // Now most of the surface conformation will be done with some
245  // residual protrusions / incursions.
246  // 2. find any segments of dual edges outside the geometry. Shoot
247  // ray from Delaunay vertex to middle of this segment and introduce
248  // point pairs. This will handle e.g.
249 
250  // protruding section of face:
251  //
252  // internal
253  // \ /
254  // -+-----------+-- boundary
255  // \ /
256  // --------
257  //
258  // Shoot ray and find intersection with outside segment (x) and
259  // introduce point pair (..)
260  //
261  // |
262  // \ . /
263  // -+-----|-----+-- boundary
264  // \ . /
265  // ---x----
266 
267  // Find vertices near boundaries to speed up subsequent checks.
268  label countNearBoundaryVertices = findVerticesNearBoundaries();
269 
270  Info<< " Vertices marked as being near a boundary: "
271  << returnReduce(countNearBoundaryVertices, sumOp<label>())
272  << " (estimated)" << endl;
273 
274  timeCheck("After set near boundary");
275 
276  const scalar edgeSearchDistCoeffSqr =
277  foamyHexMeshControls().edgeSearchDistCoeffSqr();
278 
279  const scalar surfacePtReplaceDistCoeffSqr =
280  foamyHexMeshControls().surfacePtReplaceDistCoeffSqr();
281 
282  const label AtoV = label(6/Foam::pow(scalar(number_of_vertices()), 3));
283 
284  // Initial surface protrusion conformation - nearest surface point
285  {
286  pointIndexHitAndFeatureDynList featureEdgeHits(AtoV/4);
287  pointIndexHitAndFeatureDynList surfaceHits(AtoV);
288  DynamicList<label> edgeToTreeShape(AtoV/4);
289  DynamicList<label> surfaceToTreeShape(AtoV);
290 
291  Map<scalar> surfacePtToEdgePtDist(AtoV/4);
292 
293  for
294  (
295  Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
296  vit != finite_vertices_end();
297  vit++
298  )
299  {
300  if (vit->nearBoundary())
301  {
302  pointIndexHitAndFeatureDynList surfaceIntersections(AtoV);
303 
304  if
305  (
306  dualCellSurfaceAllIntersections
307  (
308  vit,
309  surfaceIntersections
310  )
311  )
312  {
313  // meshTools::writeOBJ(Pout, vert);
314  // meshTools::writeOBJ(Pout, surfHit.hitPoint());
315  // Pout<< "l cr0 cr1" << endl;
316 
317  addSurfaceAndEdgeHits
318  (
319  topoint(vit->point()),
320  surfaceIntersections,
321  surfacePtReplaceDistCoeffSqr,
322  edgeSearchDistCoeffSqr,
323  surfaceHits,
324  featureEdgeHits,
325  surfaceToTreeShape,
326  edgeToTreeShape,
327  surfacePtToEdgePtDist,
328  true
329  );
330  }
331  else
332  {
333  vit->setInternal();
334  countNearBoundaryVertices--;
335  }
336  }
337  }
338 
339  Info<< " Vertices marked as being near a boundary: "
340  << returnReduce(countNearBoundaryVertices, sumOp<label>())
341  << " (after dual surface intersection)" << endl;
342 
343  label nVerts = number_of_vertices();
344  label nSurfHits = surfaceHits.size();
345  label nFeatEdHits = featureEdgeHits.size();
346 
347  if (Pstream::parRun())
348  {
349  reduce(nVerts, sumOp<label>());
350  reduce(nSurfHits, sumOp<label>());
351  reduce(nFeatEdHits, sumOp<label>());
352  }
353 
354  Info<< nl << "Initial conformation" << nl
355  << " Number of vertices " << nVerts << nl
356  << " Number of surface hits " << nSurfHits << nl
357  << " Number of edge hits " << nFeatEdHits
358  << endl;
359 
360  // In parallel, synchronise the surface trees
361  if (Pstream::parRun())
362  {
363  synchroniseSurfaceTrees(surfaceToTreeShape, surfaceHits);
364  }
365 
366  DynamicList<Vb> pts(2*surfaceHits.size() + 3*featureEdgeHits.size());
367 
368  insertSurfacePointPairs
369  (
370  surfaceHits,
371  "surfaceConformationLocations_initial.obj",
372  pts
373  );
374 
375  // In parallel, synchronise the edge trees
376  if (Pstream::parRun())
377  {
378  synchroniseEdgeTrees(edgeToTreeShape, featureEdgeHits);
379  }
380 
381  insertEdgePointGroups
382  (
383  featureEdgeHits,
384  "edgeConformationLocations_initial.obj",
385  pts
386  );
387 
388  pts.shrink();
389 
390  Map<label> oldToNewIndices = insertPointPairs(pts, true, true);
391 
392  // Re-index the point pairs
393  ptPairs_.reIndex(oldToNewIndices);
394 
395  //writePointPairs("pointPairs_initial.obj");
396 
397  // Remove location from surface/edge tree
398 
399  timeCheck("After initial conformation");
400 
401  initialTotalHits = nSurfHits + nFeatEdHits;
402  }
403 
404  // Remember which vertices were referred to each processor so only updates
405  // are sent.
406  PtrList<labelPairHashSet> referralVertices(Pstream::nProcs());
407 
408  // Store the vertices that have been received and added from each processor
409  // already so that there is no attempt to add them more than once.
410  autoPtr<labelPairHashSet> receivedVertices;
411 
412  if (Pstream::parRun())
413  {
414  forAll(referralVertices, proci)
415  {
416  if (proci != Pstream::myProcNo())
417  {
418  referralVertices.set
419  (
420  proci,
421  new labelPairHashSet(number_of_vertices()/Pstream::nProcs())
422  );
423  }
424  }
425 
426  receivedVertices.reset
427  (
428  new labelPairHashSet(number_of_vertices()/Pstream::nProcs())
429  );
430 
431  // Build the parallel interface the initial surface conformation
432  sync
433  (
434  decomposition_().procBounds(),
435  referralVertices,
436  receivedVertices()
437  );
438  }
439 
440  label iterationNo = 0;
441 
442  label maxIterations = foamyHexMeshControls().maxConformationIterations();
443 
444  scalar iterationToInitialHitRatioLimit =
445  foamyHexMeshControls().iterationToInitialHitRatioLimit();
446 
447  label hitLimit = label(iterationToInitialHitRatioLimit*initialTotalHits);
448 
449  Info<< nl << "Stopping iterations when: " << nl
450  << " total number of hits drops below "
451  << iterationToInitialHitRatioLimit
452  << " of initial hits (" << hitLimit << ")" << nl
453  << " or " << nl
454  << " maximum number of iterations (" << maxIterations
455  << ") is reached"
456  << endl;
457 
458  // Set totalHits to a large enough positive value to enter the while loop on
459  // the first iteration
460  label totalHits = initialTotalHits;
461 
462  while
463  (
464  totalHits > 0
465  && totalHits >= hitLimit
466  && iterationNo < maxIterations
467  )
468  {
469  pointIndexHitAndFeatureDynList surfaceHits(0.5*AtoV);
470  pointIndexHitAndFeatureDynList featureEdgeHits(0.25*AtoV);
471  DynamicList<label> surfaceToTreeShape(AtoV/2);
472  DynamicList<label> edgeToTreeShape(AtoV/4);
473 
474  Map<scalar> surfacePtToEdgePtDist;
475 
476  for
477  (
478  Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
479  vit != finite_vertices_end();
480  ++vit
481  )
482  {
483  // The initial surface conformation has already identified the
484  // nearBoundary set of vertices. Previously inserted boundary
485  // points and referred internal vertices from other processors can
486  // also generate protrusions and must be assessed too.
487  if
488  (
489  vit->nearBoundary()
490  || vit->internalBoundaryPoint()
491  || (vit->internalOrBoundaryPoint() && vit->referred())
492  )
493  {
494  pointIndexHitAndFeatureDynList surfaceIntersections(0.5*AtoV);
495 
496  pointIndexHit surfHit;
497  label hitSurface;
498 
499  // Find segments of dual face outside the geometry and find the
500  // the middle of this
501  dualCellLargestSurfaceProtrusion(vit, surfHit, hitSurface);
502 
503  if (surfHit.hit())
504  {
505  surfaceIntersections.append
506  (
507  pointIndexHitAndFeature(surfHit, hitSurface)
508  );
509 
510  addSurfaceAndEdgeHits
511  (
512  topoint(vit->point()),
513  surfaceIntersections,
514  surfacePtReplaceDistCoeffSqr,
515  edgeSearchDistCoeffSqr,
516  surfaceHits,
517  featureEdgeHits,
518  surfaceToTreeShape,
519  edgeToTreeShape,
520  surfacePtToEdgePtDist,
521  false
522  );
523  }
524  else
525  {
526  // No surface hit detected so if internal then don't check
527  // again
528  if (vit->nearBoundary())
529  {
530  vit->setInternal();
531  }
532  }
533  }
534  else if
535  (
536  vit->externalBoundaryPoint()
537  || (vit->externalBoundaryPoint() && vit->referred())
538  )
539  {
540  pointIndexHitAndFeatureDynList surfaceIntersections(0.5*AtoV);
541 
542  pointIndexHit surfHit;
543  label hitSurface;
544 
545  // Detect slave (external vertices) whose dual face incurs
546  // into nearby (other than originating) geometry
547  dualCellLargestSurfaceIncursion(vit, surfHit, hitSurface);
548 
549  if (surfHit.hit())
550  {
551  surfaceIntersections.append
552  (
553  pointIndexHitAndFeature(surfHit, hitSurface)
554  );
555 
556  addSurfaceAndEdgeHits
557  (
558  topoint(vit->point()),
559  surfaceIntersections,
560  surfacePtReplaceDistCoeffSqr,
561  edgeSearchDistCoeffSqr,
562  surfaceHits,
563  featureEdgeHits,
564  surfaceToTreeShape,
565  edgeToTreeShape,
566  surfacePtToEdgePtDist,
567  false
568  );
569  }
570  }
571  }
572 
573  label nVerts = number_of_vertices();
574  label nSurfHits = surfaceHits.size();
575  label nFeatEdHits = featureEdgeHits.size();
576 
577  if (Pstream::parRun())
578  {
579  reduce(nVerts, sumOp<label>());
580  reduce(nSurfHits, sumOp<label>());
581  reduce(nFeatEdHits, sumOp<label>());
582  }
583 
584  Info<< nl << "Conformation iteration " << iterationNo << nl
585  << " Number of vertices " << nVerts << nl
586  << " Number of surface hits " << nSurfHits << nl
587  << " Number of edge hits " << nFeatEdHits
588  << endl;
589 
590  totalHits = nSurfHits + nFeatEdHits;
591 
592  label nNotInserted = 0;
593 
594  if (totalHits > 0)
595  {
596  // In parallel, synchronise the surface trees
597  if (Pstream::parRun())
598  {
599  nNotInserted +=
600  synchroniseSurfaceTrees(surfaceToTreeShape, surfaceHits);
601  }
602 
603  DynamicList<Vb> pts
604  (
605  2*surfaceHits.size() + 3*featureEdgeHits.size()
606  );
607 
608  insertSurfacePointPairs
609  (
610  surfaceHits,
611  "surfaceConformationLocations_" + name(iterationNo) + ".obj",
612  pts
613  );
614 
615  // In parallel, synchronise the edge trees
616  if (Pstream::parRun())
617  {
618  nNotInserted +=
619  synchroniseEdgeTrees(edgeToTreeShape, featureEdgeHits);
620  }
621 
622  insertEdgePointGroups
623  (
624  featureEdgeHits,
625  "edgeConformationLocations_" + name(iterationNo) + ".obj",
626  pts
627  );
628 
629  pts.shrink();
630 
631  Map<label> oldToNewIndices = insertPointPairs(pts, true, true);
632 
633  // Reindex the point pairs
634  ptPairs_.reIndex(oldToNewIndices);
635 
636  //writePointPairs("pointPairs_" + name(iterationNo) + ".obj");
637 
638  if (Pstream::parRun())
639  {
640  sync
641  (
642  decomposition_().procBounds(),
643  referralVertices,
644  receivedVertices()
645  );
646  }
647  }
648 
649  timeCheck("Conformation iteration " + name(iterationNo));
650 
651  iterationNo++;
652 
653  if (iterationNo == maxIterations)
654  {
656  << "Maximum surface conformation iterations ("
657  << maxIterations << ") reached." << endl;
658  }
659 
660  if (totalHits <= nNotInserted)
661  {
662  Info<< nl << "Total hits (" << totalHits
663  << ") less than number of failed insertions (" << nNotInserted
664  << "), stopping iterations" << endl;
665  break;
666  }
667 
668  if (totalHits < hitLimit)
669  {
670  Info<< nl << "Total hits (" << totalHits
671  << ") less than limit (" << hitLimit
672  << "), stopping iterations" << endl;
673  }
674  }
675 
676  edgeLocationTreePtr_.clear();
677  surfacePtLocationTreePtr_.clear();
678 }
679 
680 
681 Foam::label Foam::conformalVoronoiMesh::synchroniseSurfaceTrees
682 (
683  const DynamicList<label>& surfaceToTreeShape,
684  pointIndexHitAndFeatureList& surfaceHits
685 )
686 {
687  Info<< " Surface tree synchronisation" << endl;
688 
689  pointIndexHitAndFeatureDynList synchronisedSurfLocations
690  (
691  surfaceHits.size()
692  );
693 
694  List<pointIndexHitAndFeatureDynList> procSurfLocations(Pstream::nProcs());
695 
696  procSurfLocations[Pstream::myProcNo()] = surfaceHits;
697 
698  Pstream::gatherList(procSurfLocations);
699  Pstream::scatterList(procSurfLocations);
700 
701  List<labelHashSet> hits(Pstream::nProcs());
702 
703  label nStoppedInsertion = 0;
704 
705  // Do the nearness tests here
706  for (const int proci : Pstream::allProcs())
707  {
708  // Skip own points
709  if (proci >= Pstream::myProcNo())
710  {
711  continue;
712  }
713 
714  const pointIndexHitAndFeatureList& otherSurfEdges =
715  procSurfLocations[proci];
716 
717  forAll(otherSurfEdges, peI)
718  {
719  const Foam::point& pt = otherSurfEdges[peI].first().hitPoint();
720 
721  pointIndexHit nearest;
722  pointIsNearSurfaceLocation(pt, nearest);
723 
724  pointIndexHit nearestEdge;
725  pointIsNearFeatureEdgeLocation(pt, nearestEdge);
726 
727  if (nearest.hit() || nearestEdge.hit())
728  {
729  nStoppedInsertion++;
730  hits[proci].insert(peI);
731  }
732  }
733  }
734 
735  Pstream::listCombineGather(hits, plusEqOp<labelHashSet>());
737 
738  forAll(surfaceHits, eI)
739  {
740  if (!hits[Pstream::myProcNo()].found(eI))
741  {
742  synchronisedSurfLocations.append(surfaceHits[eI]);
743  }
744  else
745  {
746  surfacePtLocationTreePtr_().remove(surfaceToTreeShape[eI]);
747  }
748  }
749 
750 // forAll(synchronisedSurfLocations, pI)
751 // {
752 // appendToSurfacePtTree
753 // (
754 // synchronisedSurfLocations[pI].first().hitPoint()
755 // );
756 // }
757 
758  const label nNotInserted = returnReduce(nStoppedInsertion, sumOp<label>());
759 
760  Info<< " Not inserting total of " << nNotInserted << " locations"
761  << endl;
762 
763  surfaceHits = synchronisedSurfLocations;
764 
765  return nNotInserted;
766 }
767 
768 
769 Foam::label Foam::conformalVoronoiMesh::synchroniseEdgeTrees
770 (
771  const DynamicList<label>& edgeToTreeShape,
772  pointIndexHitAndFeatureList& featureEdgeHits
773 )
774 {
775  Info<< " Edge tree synchronisation" << endl;
776 
777  pointIndexHitAndFeatureDynList synchronisedEdgeLocations
778  (
779  featureEdgeHits.size()
780  );
781 
782  List<pointIndexHitAndFeatureDynList> procEdgeLocations(Pstream::nProcs());
783 
784  procEdgeLocations[Pstream::myProcNo()] = featureEdgeHits;
785 
786  Pstream::gatherList(procEdgeLocations);
787  Pstream::scatterList(procEdgeLocations);
788 
789  List<labelHashSet> hits(Pstream::nProcs());
790 
791  label nStoppedInsertion = 0;
792 
793  // Do the nearness tests here
794  for (const int proci : Pstream::allProcs())
795  {
796  // Skip own points
797  if (proci >= Pstream::myProcNo())
798  {
799  continue;
800  }
801 
802  pointIndexHitAndFeatureList& otherProcEdges = procEdgeLocations[proci];
803 
804  forAll(otherProcEdges, peI)
805  {
806  const Foam::point& pt = otherProcEdges[peI].first().hitPoint();
807 
808  pointIndexHit nearest;
809  pointIsNearFeatureEdgeLocation(pt, nearest);
810 
811  if (nearest.hit())
812  {
813 // Pout<< "Not inserting " << peI << " " << pt << " "
814 // << nearest.rawPoint() << " on proc " << proci
815 // << ", near edge = " << nearest
816 // << " near ftPt = "<< info
817 // << " " << featureEdgeExclusionDistanceSqr(pt)
818 // << endl;
819 
820  nStoppedInsertion++;
821  hits[proci].insert(peI);
822  }
823  }
824  }
825 
826  Pstream::listCombineGather(hits, plusEqOp<labelHashSet>());
828 
829  forAll(featureEdgeHits, eI)
830  {
831  if (!hits[Pstream::myProcNo()].found(eI))
832  {
833  synchronisedEdgeLocations.append(featureEdgeHits[eI]);
834  }
835  else
836  {
837  edgeLocationTreePtr_().remove(edgeToTreeShape[eI]);
838  }
839  }
840 
841 // forAll(synchronisedEdgeLocations, pI)
842 // {
843 // appendToEdgeLocationTree
844 // (
845 // synchronisedEdgeLocations[pI].first().hitPoint()
846 // );
847 // }
848 
849  const label nNotInserted = returnReduce(nStoppedInsertion, sumOp<label>());
850 
851  Info<< " Not inserting total of " << nNotInserted << " locations"
852  << endl;
853 
854  featureEdgeHits = synchronisedEdgeLocations;
855 
856  return nNotInserted;
857 }
858 
859 
860 bool Foam::conformalVoronoiMesh::surfaceLocationConformsToInside
861 (
862  const pointIndexHitAndFeature& info
863 ) const
864 {
865  if (info.first().hit())
866  {
867  vectorField norm(1);
868 
869  geometryToConformTo_.getNormal
870  (
871  info.second(),
872  List<pointIndexHit>(1, info.first()),
873  norm
874  );
875 
876  const vector& n = norm[0];
877 
878  const scalar ppDist = pointPairDistance(info.first().hitPoint());
879 
880  const Foam::point innerPoint = info.first().hitPoint() - ppDist*n;
881 
882  if (!geometryToConformTo_.inside(innerPoint))
883  {
884  return false;
885  }
886 
887  return true;
888  }
889 
890  return false;
891 }
892 
893 
894 bool Foam::conformalVoronoiMesh::dualCellSurfaceAnyIntersection
895 (
896  const Delaunay::Finite_vertices_iterator& vit
897 ) const
898 {
899  std::list<Facet> facets;
900  incident_facets(vit, std::back_inserter(facets));
901 
902  for
903  (
904  std::list<Facet>::iterator fit=facets.begin();
905  fit != facets.end();
906  ++fit
907  )
908  {
909  if
910  (
911  is_infinite(fit->first)
912  || is_infinite(fit->first->neighbor(fit->second))
913  || !fit->first->hasInternalPoint()
914  || !fit->first->neighbor(fit->second)->hasInternalPoint()
915  )
916  {
917  continue;
918  }
919 
920  Foam::point dE0 = fit->first->dual();
921  Foam::point dE1 = fit->first->neighbor(fit->second)->dual();
922 
923  if (Pstream::parRun())
924  {
925  Foam::point& a = dE0;
926  Foam::point& b = dE1;
927 
928  bool inProc = clipLineToProc(topoint(vit->point()), a, b);
929 
930  // Check for the edge passing through a surface
931  if
932  (
933  inProc
934  && geometryToConformTo_.findSurfaceAnyIntersection(a, b)
935  )
936  {
937  return true;
938  }
939  }
940  else
941  {
942  if (geometryToConformTo_.findSurfaceAnyIntersection(dE0, dE1))
943  {
944  return true;
945  }
946  }
947  }
948 
949  return false;
950 }
951 
952 
953 bool Foam::conformalVoronoiMesh::dualCellSurfaceAllIntersections
954 (
955  const Delaunay::Finite_vertices_iterator& vit,
956  pointIndexHitAndFeatureDynList& infoList
957 ) const
958 {
959  bool flagIntersection = false;
960 
961  std::list<Facet> facets;
962  incident_facets(vit, std::back_inserter(facets));
963 
964  for
965  (
966  std::list<Facet>::iterator fit = facets.begin();
967  fit != facets.end();
968  ++fit
969  )
970  {
971  if
972  (
973  is_infinite(fit->first)
974  || is_infinite(fit->first->neighbor(fit->second))
975  || !fit->first->hasInternalPoint()
976  || !fit->first->neighbor(fit->second)->hasInternalPoint()
977  )
978  {
979  continue;
980  }
981 
982  // Construct the dual edge and search for intersections of the edge
983  // with the surface
984  Foam::point dE0 = fit->first->dual();
985  Foam::point dE1 = fit->first->neighbor(fit->second)->dual();
986 
987  pointIndexHit infoIntersection;
988  label hitSurfaceIntersection = -1;
989 
990  if (Pstream::parRun())
991  {
992  bool inProc = clipLineToProc(topoint(vit->point()), dE0, dE1);
993 
994  if (!inProc)
995  {
996  continue;
997  }
998  }
999 
1000  geometryToConformTo_.findSurfaceNearestIntersection
1001  (
1002  dE0,
1003  dE1,
1004  infoIntersection,
1005  hitSurfaceIntersection
1006  );
1007 
1008  if (infoIntersection.hit())
1009  {
1010  vectorField norm(1);
1011 
1012  geometryToConformTo_.getNormal
1013  (
1014  hitSurfaceIntersection,
1015  List<pointIndexHit>(1, infoIntersection),
1016  norm
1017  );
1018 
1019  const vector& n = norm[0];
1020 
1021  pointFromPoint vertex = topoint(vit->point());
1022 
1023  const plane p(infoIntersection.hitPoint(), n);
1024 
1025  const plane::ray r(vertex, n);
1026 
1027  const scalar d = p.normalIntersect(r);
1028 
1029  Foam::point newPoint = vertex + d*n;
1030 
1031  pointIndexHitAndFeature info;
1032  geometryToConformTo_.findSurfaceNearest
1033  (
1034  newPoint,
1035  4.0*magSqr(newPoint - vertex),
1036  info.first(),
1037  info.second()
1038  );
1039 
1040  bool rejectPoint = false;
1041 
1042  if (!surfaceLocationConformsToInside(info))
1043  {
1044  rejectPoint = true;
1045  }
1046 
1047  if (!rejectPoint && info.first().hit())
1048  {
1049  if (!infoList.empty())
1050  {
1051  forAll(infoList, hitI)
1052  {
1053  // Reject point if the point is already added
1054  if
1055  (
1056  infoList[hitI].first().index()
1057  == info.first().index()
1058  )
1059  {
1060  rejectPoint = true;
1061  break;
1062  }
1063 
1064  const Foam::point& p
1065  = infoList[hitI].first().hitPoint();
1066 
1067  const scalar separationDistance =
1068  mag(p - info.first().hitPoint());
1069 
1070  const scalar minSepDist =
1071  sqr
1072  (
1073  foamyHexMeshControls().removalDistCoeff()
1074  *targetCellSize(p)
1075  );
1076 
1077  // Reject the point if it is too close to another
1078  // surface point.
1079  // Could merge the points?
1080  if (separationDistance < minSepDist)
1081  {
1082  rejectPoint = true;
1083  break;
1084  }
1085  }
1086  }
1087  }
1088 
1089  // The normal ray from the vertex will not always result in a hit
1090  // because another surface may be in the way.
1091  if (!rejectPoint && info.first().hit())
1092  {
1093  flagIntersection = true;
1094  infoList.append(info);
1095  }
1096  }
1097  }
1098 
1099  return flagIntersection;
1100 }
1101 
1102 
1103 bool Foam::conformalVoronoiMesh::clipLineToProc
1104 (
1105  const Foam::point& pt,
1106  Foam::point& a,
1107  Foam::point& b
1108 ) const
1109 {
1110  bool inProc = false;
1111 
1112  pointIndexHit findAnyIntersection = decomposition_().findLine(a, b);
1113 
1114  if (!findAnyIntersection.hit())
1115  {
1116  pointIndexHit info = decomposition_().findLine(a, pt);
1117 
1118  if (!info.hit())
1119  {
1120  inProc = true;
1121  }
1122  else
1123  {
1124  inProc = false;
1125  }
1126  }
1127  else
1128  {
1129  pointIndexHit info = decomposition_().findLine(a, pt);
1130 
1131  if (!info.hit())
1132  {
1133  inProc = true;
1134  b = findAnyIntersection.hitPoint();
1135  }
1136  else
1137  {
1138  inProc = true;
1139  a = findAnyIntersection.hitPoint();
1140  }
1141  }
1142 
1143  return inProc;
1144 }
1145 
1146 
1147 void Foam::conformalVoronoiMesh::dualCellLargestSurfaceProtrusion
1148 (
1149  const Delaunay::Finite_vertices_iterator& vit,
1150  pointIndexHit& surfHitLargest,
1151  label& hitSurfaceLargest
1152 ) const
1153 {
1154  // Set no-hit data
1155  surfHitLargest = pointIndexHit();
1156  hitSurfaceLargest = -1;
1157 
1158  std::list<Facet> facets;
1159  finite_incident_facets(vit, std::back_inserter(facets));
1160 
1161  pointFromPoint vert = topoint(vit->point());
1162 
1163  scalar maxProtrusionDistance = maxSurfaceProtrusion(vert);
1164 
1165  for
1166  (
1167  std::list<Facet>::iterator fit = facets.begin();
1168  fit != facets.end();
1169  ++fit
1170  )
1171  {
1172  Cell_handle c1 = fit->first;
1173  Cell_handle c2 = fit->first->neighbor(fit->second);
1174 
1175  if
1176  (
1177  is_infinite(c1) || is_infinite(c2)
1178  || (
1179  !c1->internalOrBoundaryDualVertex()
1180  || !c2->internalOrBoundaryDualVertex()
1181  )
1182  || !c1->real() || !c2->real()
1183  )
1184  {
1185  continue;
1186  }
1187 
1188 // Foam::point endPt = 0.5*(c1->dual() + c2->dual());
1189  Foam::point endPt = c1->dual();
1190 
1191  if (magSqr(vert - c1->dual()) < magSqr(vert - c2->dual()))
1192  {
1193  endPt = c2->dual();
1194  }
1195 
1196  if
1197  (
1198  magSqr(vert - endPt)
1199  > magSqr(geometryToConformTo().globalBounds().mag())
1200  )
1201  {
1202  continue;
1203  }
1204 
1205  pointIndexHit surfHit;
1206  label hitSurface;
1207 
1208  geometryToConformTo_.findSurfaceNearestIntersection
1209  (
1210  vert,
1211  endPt,
1212  surfHit,
1213  hitSurface
1214  );
1215 
1216  if (surfHit.hit())
1217  {
1218  vectorField norm(1);
1219 
1220  allGeometry_[hitSurface].getNormal
1221  (
1222  List<pointIndexHit>(1, surfHit),
1223  norm
1224  );
1225 
1226  const vector& n = norm[0];
1227 
1228  const scalar normalProtrusionDistance
1229  (
1230  (endPt - surfHit.hitPoint()) & n
1231  );
1232 
1233  if (normalProtrusionDistance > maxProtrusionDistance)
1234  {
1235  const plane p(surfHit.hitPoint(), n);
1236 
1237  const plane::ray r(endPt, -n);
1238 
1239  const scalar d = p.normalIntersect(r);
1240 
1241  Foam::point newPoint = endPt - d*n;
1242 
1243  pointIndexHitAndFeature info;
1244  geometryToConformTo_.findSurfaceNearest
1245  (
1246  newPoint,
1247  4.0*magSqr(newPoint - endPt),
1248  info.first(),
1249  info.second()
1250  );
1251 
1252  if (info.first().hit())
1253  {
1254  if
1255  (
1256  surfaceLocationConformsToInside
1257  (
1258  pointIndexHitAndFeature(info.first(), info.second())
1259  )
1260  )
1261  {
1262  surfHitLargest = info.first();
1263  hitSurfaceLargest = info.second();
1264 
1265  maxProtrusionDistance = normalProtrusionDistance;
1266  }
1267  }
1268  }
1269  }
1270  }
1271 
1272  // Relying on short-circuit evaluation to not call for hitPoint when this
1273  // is a miss
1274  if
1275  (
1276  surfHitLargest.hit()
1277  && (
1278  Pstream::parRun()
1279  && !decomposition().positionOnThisProcessor(surfHitLargest.hitPoint())
1280  )
1281  )
1282  {
1283  // A protrusion was identified, but not penetrating on this processor,
1284  // so set no-hit data and allow the other that should have this point
1285  // referred to generate it.
1286  surfHitLargest = pointIndexHit();
1287  hitSurfaceLargest = -1;
1288  }
1289 }
1290 
1291 
1292 void Foam::conformalVoronoiMesh::dualCellLargestSurfaceIncursion
1293 (
1294  const Delaunay::Finite_vertices_iterator& vit,
1295  pointIndexHit& surfHitLargest,
1296  label& hitSurfaceLargest
1297 ) const
1298 {
1299  // Set no-hit data
1300  surfHitLargest = pointIndexHit();
1301  hitSurfaceLargest = -1;
1302 
1303  std::list<Facet> facets;
1304  finite_incident_facets(vit, std::back_inserter(facets));
1305 
1306  pointFromPoint vert = topoint(vit->point());
1307 
1308  scalar minIncursionDistance = -maxSurfaceProtrusion(vert);
1309 
1310  for
1311  (
1312  std::list<Facet>::iterator fit = facets.begin();
1313  fit != facets.end();
1314  ++fit
1315  )
1316  {
1317  Cell_handle c1 = fit->first;
1318  Cell_handle c2 = fit->first->neighbor(fit->second);
1319 
1320  if
1321  (
1322  is_infinite(c1) || is_infinite(c2)
1323  || (
1324  !c1->internalOrBoundaryDualVertex()
1325  || !c2->internalOrBoundaryDualVertex()
1326  )
1327  || !c1->real() || !c2->real()
1328  )
1329  {
1330  continue;
1331  }
1332 
1333 // Foam::point endPt = 0.5*(c1->dual() + c2->dual());
1334  Foam::point endPt = c1->dual();
1335 
1336  if (magSqr(vert - c1->dual()) < magSqr(vert - c2->dual()))
1337  {
1338  endPt = c2->dual();
1339  }
1340 
1341  if
1342  (
1343  magSqr(vert - endPt)
1344  > magSqr(geometryToConformTo().globalBounds().mag())
1345  )
1346  {
1347  continue;
1348  }
1349 
1350  pointIndexHit surfHit;
1351  label hitSurface;
1352 
1353  geometryToConformTo_.findSurfaceNearestIntersection
1354  (
1355  vert,
1356  endPt,
1357  surfHit,
1358  hitSurface
1359  );
1360 
1361  if (surfHit.hit())
1362  {
1363  vectorField norm(1);
1364 
1365  allGeometry_[hitSurface].getNormal
1366  (
1367  List<pointIndexHit>(1, surfHit),
1368  norm
1369  );
1370 
1371  const vector& n = norm[0];
1372 
1373  scalar normalIncursionDistance
1374  (
1375  (endPt - surfHit.hitPoint()) & n
1376  );
1377 
1378  if (normalIncursionDistance < minIncursionDistance)
1379  {
1380  const plane p(surfHit.hitPoint(), n);
1381 
1382  const plane::ray r(endPt, n);
1383 
1384  const scalar d = p.normalIntersect(r);
1385 
1386  Foam::point newPoint = endPt + d*n;
1387 
1388  pointIndexHitAndFeature info;
1389  geometryToConformTo_.findSurfaceNearest
1390  (
1391  newPoint,
1392  4.0*magSqr(newPoint - endPt),
1393  info.first(),
1394  info.second()
1395  );
1396 
1397  if (info.first().hit())
1398  {
1399  if
1400  (
1401  surfaceLocationConformsToInside
1402  (
1403  pointIndexHitAndFeature(info.first(), info.second())
1404  )
1405  )
1406  {
1407  surfHitLargest = info.first();
1408  hitSurfaceLargest = info.second();
1409 
1410  minIncursionDistance = normalIncursionDistance;
1411  }
1412  }
1413  }
1414  }
1415  }
1416 
1417  // Relying on short-circuit evaluation to not call for hitPoint when this
1418  // is a miss
1419  if
1420  (
1421  surfHitLargest.hit()
1422  && (
1423  Pstream::parRun()
1424  && !decomposition().positionOnThisProcessor(surfHitLargest.hitPoint())
1425  )
1426  )
1427  {
1428  // A protrusion was identified, but not penetrating on this processor,
1429  // so set no-hit data and allow the other that should have this point
1430  // referred to generate it.
1431  surfHitLargest = pointIndexHit();
1432  hitSurfaceLargest = -1;
1433  }
1434 }
1435 
1436 
1437 void Foam::conformalVoronoiMesh::reportProcessorOccupancy()
1438 {
1439  for
1440  (
1441  Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1442  vit != finite_vertices_end();
1443  vit++
1444  )
1445  {
1446  if (vit->real())
1447  {
1448  if
1449  (
1450  Pstream::parRun()
1451  && !decomposition().positionOnThisProcessor(topoint(vit->point()))
1452  )
1453  {
1454  Pout<< topoint(vit->point()) << " is not on this processor "
1455  << endl;
1456  }
1457  }
1458  }
1459 }
1460 
1461 
1462 //void Foam::conformalVoronoiMesh::reportSurfaceConformationQuality()
1463 //{
1464 // Info<< nl << "Check surface conformation quality" << endl;
1465 //
1466 // for
1467 // (
1468 // Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1469 // vit != finite_vertices_end();
1470 // vit++
1471 // )
1472 // {
1473 // if (vit->internalOrBoundaryPoint())
1474 // {
1475 // Foam::point vert(topoint(vit->point()));
1476 // pointIndexHit surfHit;
1477 // label hitSurface;
1478 //
1479 // dualCellLargestSurfaceProtrusion(vit, surfHit, hitSurface);
1480 //
1481 // if (surfHit.hit())
1482 // {
1483 // Pout<< nl << "Residual penetration: " << nl
1484 // << vit->index() << nl
1485 // << vit->type() << nl
1486 // << vit->ppMaster() << nl
1487 // << "nearFeaturePt "
1488 // << nearFeaturePt(surfHit.hitPoint()) << nl
1489 // << vert << nl
1490 // << surfHit.hitPoint()
1491 // << endl;
1492 // }
1493 // }
1494 // }
1495 //
1496 // {
1497 // // Assess close surface points
1498 //
1499 // setVertexSizeAndAlignment();
1500 //
1501 // for
1502 // (
1503 // Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1504 // vit != finite_vertices_end();
1505 // vit++
1506 // )
1507 // {
1508 // if (vit->ppMaster())
1509 // {
1510 // std::list<Vertex_handle> adjacentVertices;
1511 //
1512 // adjacent_vertices(vit, std::back_inserter(adjacentVertices));
1513 //
1514 // Foam::point pt = topoint(vit->point());
1515 //
1516 // // Pout<< nl << "vit: " << vit->index() << " "
1517 // // << topoint(vit->point())
1518 // // << endl;
1519 //
1520 // // Pout<< adjacentVertices.size() << endl;
1521 //
1522 // for
1523 // (
1524 // std::list<Vertex_handle>::iterator
1525 // avit = adjacentVertices.begin();
1526 // avit != adjacentVertices.end();
1527 // ++avit
1528 // )
1529 // {
1530 // Vertex_handle avh = *avit;
1531 //
1532 // // The lower indexed vertex will perform the assessment
1533 // if
1534 // (
1535 // avh->ppMaster()
1536 // && vit->index() < avh->index()
1537 // && vit->type() != avh->type()
1538 // )
1539 // {
1540 // scalar targetSize = 0.2*averageAnyCellSize(vit, avh);
1541 //
1542 // // Pout<< "diff " << mag(pt - topoint(avh->point()))
1543 // // << " " << targetSize << endl;
1544 //
1545 // if
1546 // (
1547 // magSqr(pt - topoint(avh->point()))
1548 // < sqr(targetSize)
1549 // )
1550 // {
1551 // Pout<< nl << "vit: " << vit->index() << " "
1552 // << topoint(vit->point())
1553 // << endl;
1554 //
1555 // Pout<< " adjacent too close: "
1556 // << avh->index() << " "
1557 // << topoint(avh->point())
1558 // << endl;
1559 // }
1560 // }
1561 // }
1562 // }
1563 // }
1564 // }
1565 //}
1566 
1567 void Foam::conformalVoronoiMesh::limitDisplacement
1568 (
1569  const Delaunay::Finite_vertices_iterator& vit,
1570  vector& displacement,
1571  label callCount
1572 ) const
1573 {
1574  callCount++;
1575 
1576  // Do not allow infinite recursion
1577  if (callCount > 7)
1578  {
1579  displacement = Zero;
1580  return;
1581  }
1582 
1583  pointFromPoint pt = topoint(vit->point());
1584  Foam::point dispPt = pt + displacement;
1585 
1586  bool limit = false;
1587 
1588  pointIndexHit surfHit;
1589  label hitSurface;
1590 
1591  if (!geometryToConformTo_.globalBounds().contains(dispPt))
1592  {
1593  // If dispPt is outside bounding box then displacement cuts boundary
1594  limit = true;
1595  }
1596  else if (geometryToConformTo_.findSurfaceAnyIntersection(pt, dispPt))
1597  {
1598  // Full surface penetration test
1599  limit = true;
1600  }
1601  else
1602  {
1603  // Testing if the displaced position is too close to the surface.
1604  // Within twice the local surface point pair insertion distance is
1605  // considered "too close"
1606 
1607  scalar searchDistanceSqr = sqr
1608  (
1609  2*vit->targetCellSize()
1610  *foamyHexMeshControls().pointPairDistanceCoeff()
1611  );
1612 
1613  geometryToConformTo_.findSurfaceNearest
1614  (
1615  dispPt,
1616  searchDistanceSqr,
1617  surfHit,
1618  hitSurface
1619  );
1620 
1621  if (surfHit.hit())
1622  {
1623  limit = true;
1624 
1625  if (magSqr(pt - surfHit.hitPoint()) <= searchDistanceSqr)
1626  {
1627  // Cannot limit displacement, point closer than tolerance
1628  displacement = Zero;
1629  return;
1630  }
1631  }
1632  }
1633 
1634  if (limit)
1635  {
1636  // Halve the displacement and call this function again. Will continue
1637  // recursively until the displacement is small enough.
1638 
1639  displacement *= 0.5;
1640 
1641  limitDisplacement(vit, displacement, callCount);
1642  }
1643 }
1644 
1645 
1646 Foam::scalar Foam::conformalVoronoiMesh::angleBetweenSurfacePoints
1647 (
1648  Foam::point pA,
1649  Foam::point pB
1650 ) const
1651 {
1652  pointIndexHit pAhit;
1653  label pAsurfaceHit = -1;
1654 
1655  const scalar searchDist = 5.0*targetCellSize(pA);
1656 
1657  geometryToConformTo_.findSurfaceNearest
1658  (
1659  pA,
1660  searchDist,
1661  pAhit,
1662  pAsurfaceHit
1663  );
1664 
1665  if (!pAhit.hit())
1666  {
1668  }
1669 
1670  vectorField norm(1);
1671 
1672  allGeometry_[pAsurfaceHit].getNormal
1673  (
1674  List<pointIndexHit>(1, pAhit),
1675  norm
1676  );
1677 
1678  const vector nA = norm[0];
1679 
1680  pointIndexHit pBhit;
1681  label pBsurfaceHit = -1;
1682 
1683  geometryToConformTo_.findSurfaceNearest
1684  (
1685  pB,
1686  searchDist,
1687  pBhit,
1688  pBsurfaceHit
1689  );
1690 
1691  if (!pBhit.hit())
1692  {
1694  }
1695 
1696  allGeometry_[pBsurfaceHit].getNormal
1697  (
1698  List<pointIndexHit>(1, pBhit),
1699  norm
1700  );
1701 
1702  const vector nB = norm[0];
1703 
1704  return vectorTools::cosPhi(nA, nB);
1705 }
1706 
1707 
1708 bool Foam::conformalVoronoiMesh::nearSurfacePoint
1709 (
1710  pointIndexHitAndFeature& pHit
1711 ) const
1712 {
1713  const Foam::point& pt = pHit.first().hitPoint();
1714 
1715  pointIndexHit closePoint;
1716  const bool closeToSurfacePt = pointIsNearSurfaceLocation(pt, closePoint);
1717 
1718  if
1719  (
1720  closeToSurfacePt
1721  && (
1722  magSqr(pt - closePoint.hitPoint())
1723  > sqr(pointPairDistance(pt))
1724  )
1725  )
1726  {
1727  const scalar cosAngle =
1728  angleBetweenSurfacePoints(pt, closePoint.hitPoint());
1729 
1730  // TODO: make this tolerance run-time selectable?
1731  if (cosAngle < searchAngleOppositeSurface)
1732  {
1733  pointIndexHit pCloseHit;
1734  label pCloseSurfaceHit = -1;
1735 
1736  const scalar searchDist = targetCellSize(closePoint.hitPoint());
1737 
1738  geometryToConformTo_.findSurfaceNearest
1739  (
1740  closePoint.hitPoint(),
1741  searchDist,
1742  pCloseHit,
1743  pCloseSurfaceHit
1744  );
1745 
1746  vectorField norm(1);
1747 
1748  allGeometry_[pCloseSurfaceHit].getNormal
1749  (
1750  List<pointIndexHit>(1, pCloseHit),
1751  norm
1752  );
1753 
1754  const vector& nA = norm[0];
1755 
1756  pointIndexHit oppositeHit;
1757  label oppositeSurfaceHit = -1;
1758 
1759  geometryToConformTo_.findSurfaceNearestIntersection
1760  (
1761  closePoint.hitPoint() + 0.5*pointPairDistance(pt)*nA,
1762  closePoint.hitPoint() + 5*targetCellSize(pt)*nA,
1763  oppositeHit,
1764  oppositeSurfaceHit
1765  );
1766 
1767  if (oppositeHit.hit())
1768  {
1769  // Replace point
1770  pHit.first() = oppositeHit;
1771  pHit.second() = oppositeSurfaceHit;
1772 
1773  return !closeToSurfacePt;
1774  }
1775  }
1776  }
1777 
1778  return closeToSurfacePt;
1779 }
1780 
1781 
1782 bool Foam::conformalVoronoiMesh::appendToSurfacePtTree
1783 (
1784  const Foam::point& pt
1785 ) const
1786 {
1787  label startIndex = existingSurfacePtLocations_.size();
1788 
1789  existingSurfacePtLocations_.append(pt);
1790 
1791  label endIndex = existingSurfacePtLocations_.size();
1792 
1793  return surfacePtLocationTreePtr_().insert(startIndex, endIndex);
1794 }
1795 
1796 
1797 bool Foam::conformalVoronoiMesh::appendToEdgeLocationTree
1798 (
1799  const Foam::point& pt
1800 ) const
1801 {
1802  label startIndex = existingEdgeLocations_.size();
1803 
1804  existingEdgeLocations_.append(pt);
1805 
1806  label endIndex = existingEdgeLocations_.size();
1807 
1808  return edgeLocationTreePtr_().insert(startIndex, endIndex);
1809 }
1810 
1811 
1813 Foam::conformalVoronoiMesh::nearestFeatureEdgeLocations
1814 (
1815  const Foam::point& pt
1816 ) const
1817 {
1818  const scalar exclusionRangeSqr = featureEdgeExclusionDistanceSqr(pt);
1819 
1820  labelList elems
1821  = edgeLocationTreePtr_().findSphere(pt, exclusionRangeSqr);
1822 
1823  DynamicList<pointIndexHit> dynPointHit;
1824 
1825  forAll(elems, elemI)
1826  {
1827  label index = elems[elemI];
1828 
1829  const Foam::point& pointi
1830  = edgeLocationTreePtr_().shapes().shapePoints()[index];
1831 
1832  pointIndexHit nearHit(true, pointi, index);
1833 
1834  dynPointHit.append(nearHit);
1835  }
1836 
1837  return dynPointHit;
1838 }
1839 
1840 
1841 bool Foam::conformalVoronoiMesh::pointIsNearFeatureEdgeLocation
1842 (
1843  const Foam::point& pt
1844 ) const
1845 {
1846  const scalar exclusionRangeSqr = featureEdgeExclusionDistanceSqr(pt);
1847 
1848  pointIndexHit info
1849  = edgeLocationTreePtr_().findNearest(pt, exclusionRangeSqr);
1850 
1851  return info.hit();
1852 }
1853 
1854 
1855 bool Foam::conformalVoronoiMesh::pointIsNearFeatureEdgeLocation
1856 (
1857  const Foam::point& pt,
1858  pointIndexHit& info
1859 ) const
1860 {
1861  const scalar exclusionRangeSqr = featureEdgeExclusionDistanceSqr(pt);
1862 
1863  info = edgeLocationTreePtr_().findNearest(pt, exclusionRangeSqr);
1864 
1865  return info.hit();
1866 }
1867 
1868 
1869 bool Foam::conformalVoronoiMesh::pointIsNearSurfaceLocation
1870 (
1871  const Foam::point& pt
1872 ) const
1873 {
1874  pointIndexHit info;
1875 
1876  pointIsNearSurfaceLocation(pt, info);
1877 
1878  return info.hit();
1879 }
1880 
1881 
1882 bool Foam::conformalVoronoiMesh::pointIsNearSurfaceLocation
1883 (
1884  const Foam::point& pt,
1885  pointIndexHit& info
1886 ) const
1887 {
1888  const scalar exclusionRangeSqr = surfacePtExclusionDistanceSqr(pt);
1889 
1890  info = surfacePtLocationTreePtr_().findNearest(pt, exclusionRangeSqr);
1891 
1892  return info.hit();
1893 }
1894 
1895 
1896 bool Foam::conformalVoronoiMesh::nearFeatureEdgeLocation
1897 (
1898  const pointIndexHit& pHit,
1899  pointIndexHit& nearestEdgeHit
1900 ) const
1901 {
1902  const Foam::point& pt = pHit.hitPoint();
1903 
1904  const scalar exclusionRangeSqr = featureEdgeExclusionDistanceSqr(pt);
1905 
1906  bool closeToFeatureEdge =
1907  pointIsNearFeatureEdgeLocation(pt, nearestEdgeHit);
1908 
1909  if (closeToFeatureEdge)
1910  {
1911  List<pointIndexHit> nearHits = nearestFeatureEdgeLocations(pt);
1912 
1913  forAll(nearHits, elemI)
1914  {
1915  pointIndexHit& info = nearHits[elemI];
1916 
1917  // Check if the edge location that the new edge location is near to
1918  // "might" be on a different edge. If so, add it anyway.
1919  pointIndexHit edgeHit;
1920  label featureHit = -1;
1921 
1922  geometryToConformTo_.findEdgeNearest
1923  (
1924  pt,
1925  exclusionRangeSqr,
1926  edgeHit,
1927  featureHit
1928  );
1929 
1930  const extendedFeatureEdgeMesh& eMesh
1931  = geometryToConformTo_.features()[featureHit];
1932 
1933  const vector& edgeDir = eMesh.edgeDirections()[edgeHit.index()];
1934 
1935  const vector lineBetweenPoints = pt - info.hitPoint();
1936 
1937  const scalar cosAngle
1938  = vectorTools::cosPhi(edgeDir, lineBetweenPoints);
1939 
1940  // Allow the point to be added if it is almost at right angles to
1941  // the other point. Also check it is not the same point.
1942  // Info<< cosAngle<< " "
1943  // << radToDeg(acos(cosAngle)) << " "
1944  // << searchConeAngle << " "
1945  // << radToDeg(acos(searchConeAngle)) << endl;
1946 
1947  if
1948  (
1949  mag(cosAngle) < searchConeAngle
1950  && (mag(lineBetweenPoints) > pointPairDistance(pt))
1951  )
1952  {
1953  //pt = edgeHit.hitPoint();
1954  //pHit.setPoint(pt);
1955  closeToFeatureEdge = false;
1956  }
1957  else
1958  {
1959  closeToFeatureEdge = true;
1960  break;
1961  }
1962  }
1963  }
1964 
1965  return closeToFeatureEdge;
1966 }
1967 
1968 
1969 void Foam::conformalVoronoiMesh::buildEdgeLocationTree
1970 (
1971  const DynamicList<Foam::point>& existingEdgeLocations
1972 ) const
1973 {
1974  treeBoundBox overallBb
1975  (
1976  geometryToConformTo_.globalBounds().extend(rndGen_, 1e-4)
1977  );
1978 
1979  overallBb.min() -= Foam::point::uniform(ROOTVSMALL);
1980  overallBb.max() += Foam::point::uniform(ROOTVSMALL);
1981 
1982  edgeLocationTreePtr_.reset
1983  (
1984  new dynamicIndexedOctree<dynamicTreeDataPoint>
1985  (
1986  dynamicTreeDataPoint(existingEdgeLocations),
1987  overallBb, // overall search domain
1988  10, // max levels, n/a
1989  20.0, // maximum ratio of cubes v.s. cells
1990  100.0 // max. duplicity; n/a since no bounding boxes.
1991  )
1992  );
1993 }
1994 
1995 
1996 void Foam::conformalVoronoiMesh::buildSurfacePtLocationTree
1997 (
1998  const DynamicList<Foam::point>& existingSurfacePtLocations
1999 ) const
2000 {
2001  treeBoundBox overallBb
2002  (
2003  geometryToConformTo_.globalBounds().extend(rndGen_, 1e-4)
2004  );
2005 
2006  overallBb.min() -= Foam::point::uniform(ROOTVSMALL);
2007  overallBb.max() += Foam::point::uniform(ROOTVSMALL);
2008 
2009  surfacePtLocationTreePtr_.reset
2010  (
2011  new dynamicIndexedOctree<dynamicTreeDataPoint>
2012  (
2013  dynamicTreeDataPoint(existingSurfacePtLocations),
2014  overallBb, // overall search domain
2015  10, // max levels, n/a
2016  20.0, // maximum ratio of cubes v.s. cells
2017  100.0 // max. duplicity; n/a since no bounding boxes.
2018  )
2019  );
2020 }
2021 
2022 
2023 void Foam::conformalVoronoiMesh::addSurfaceAndEdgeHits
2024 (
2025  const Foam::point& vit,
2026  const pointIndexHitAndFeatureDynList& surfaceIntersections,
2027  scalar surfacePtReplaceDistCoeffSqr,
2028  scalar edgeSearchDistCoeffSqr,
2029  pointIndexHitAndFeatureDynList& surfaceHits,
2030  pointIndexHitAndFeatureDynList& featureEdgeHits,
2031  DynamicList<label>& surfaceToTreeShape,
2032  DynamicList<label>& edgeToTreeShape,
2033  Map<scalar>& surfacePtToEdgePtDist,
2034  bool firstPass
2035 ) const
2036 {
2037  const scalar cellSize = targetCellSize(vit);
2038  const scalar cellSizeSqr = sqr(cellSize);
2039 
2040  forAll(surfaceIntersections, sI)
2041  {
2042  pointIndexHitAndFeature surfHitI = surfaceIntersections[sI];
2043 
2044  bool keepSurfacePoint = true;
2045 
2046  if (!surfHitI.first().hit())
2047  {
2048  continue;
2049  }
2050 
2051  const Foam::point& surfPt = surfHitI.first().hitPoint();
2052 
2053  bool isNearFeaturePt = nearFeaturePt(surfPt);
2054 
2055  bool isNearFeatureEdge = surfacePtNearFeatureEdge(surfPt);
2056 
2057  bool isNearSurfacePt = nearSurfacePoint(surfHitI);
2058 
2059  if (isNearFeaturePt || isNearSurfacePt || isNearFeatureEdge)
2060  {
2061  keepSurfacePoint = false;
2062  }
2063 
2064  List<List<pointIndexHit>> edHitsByFeature;
2065 
2066  labelList featuresHit;
2067 
2068  const scalar searchRadiusSqr = edgeSearchDistCoeffSqr*cellSizeSqr;
2069 
2070  geometryToConformTo_.findAllNearestEdges
2071  (
2072  surfPt,
2073  searchRadiusSqr,
2074  edHitsByFeature,
2075  featuresHit
2076  );
2077 
2078  forAll(edHitsByFeature, i)
2079  {
2080  const label featureHit = featuresHit[i];
2081 
2082  List<pointIndexHit>& edHits = edHitsByFeature[i];
2083 
2084  forAll(edHits, eHitI)
2085  {
2086  pointIndexHit& edHit = edHits[eHitI];
2087 
2088  if (edHit.hit())
2089  {
2090  const Foam::point& edPt = edHit.hitPoint();
2091 
2092  if
2093  (
2094  Pstream::parRun()
2095  && !decomposition().positionOnThisProcessor(edPt)
2096  )
2097  {
2098  // Do not insert
2099  continue;
2100  }
2101 
2102  if (!nearFeaturePt(edPt))
2103  {
2104  if
2105  (
2106  magSqr(edPt - surfPt)
2107  < surfacePtReplaceDistCoeffSqr*cellSizeSqr
2108  )
2109  {
2110  // If the point is within a given distance of a
2111  // feature edge, give control to edge control points
2112  // instead, this will prevent "pits" forming.
2113 
2114  // Allow if different surfaces
2115 
2116 
2117  keepSurfacePoint = false;
2118  }
2119 
2120  pointIndexHit nearestEdgeHit;
2121 
2122  if
2123  (
2124 // !pointIsNearFeatureEdgeLocation
2125 // (
2126 // edPt,
2127 // nearestEdgeHit
2128 // )
2129  !nearFeatureEdgeLocation(edHit, nearestEdgeHit)
2130  )
2131  {
2132  appendToEdgeLocationTree(edPt);
2133 
2134  edgeToTreeShape.append
2135  (
2136  existingEdgeLocations_.size() - 1
2137  );
2138 
2139  // Do not place edge control points too close to a
2140  // feature point or existing edge control points
2141  featureEdgeHits.append
2142  (
2143  pointIndexHitAndFeature(edHit, featureHit)
2144  );
2145 
2146 // Info<< "Add " << existingEdgeLocations_.size() - 1
2147 // << " " << magSqr(edPt - surfPt) << endl;
2148 
2149  surfacePtToEdgePtDist.insert
2150  (
2151  existingEdgeLocations_.size() - 1,
2152  magSqr(edPt - surfPt)
2153  );
2154  }
2155  else if (firstPass)
2156  {
2157  label hitIndex = nearestEdgeHit.index();
2158 
2159 // Info<< "Close to " << nearestEdgeHit << endl;
2160 
2161  if
2162  (
2163  magSqr(edPt - surfPt)
2164  < surfacePtToEdgePtDist[hitIndex]
2165  )
2166  {
2167  featureEdgeHits[hitIndex] =
2168  pointIndexHitAndFeature(edHit, featureHit);
2169 
2170  existingEdgeLocations_[hitIndex] =
2171  edHit.hitPoint();
2172  surfacePtToEdgePtDist[hitIndex] =
2173  magSqr(edPt - surfPt);
2174 
2175  // Change edge location in featureEdgeHits
2176  // remove index from edge tree
2177  // reinsert new point into tree
2178  edgeLocationTreePtr_().remove(hitIndex);
2179  edgeLocationTreePtr_().insert
2180  (
2181  hitIndex,
2182  hitIndex + 1
2183  );
2184  }
2185  }
2186  }
2187  }
2188  }
2189  }
2190 
2191  if (keepSurfacePoint)
2192  {
2193  surfaceHits.append(surfHitI);
2194  appendToSurfacePtTree(surfPt);
2195  surfaceToTreeShape.append(existingSurfacePtLocations_.size() - 1);
2196 
2197 // addedPoints.write(surfPt);
2198  }
2199  else
2200  {
2201 // removedPoints.write(surfPt);
2202  }
2203  }
2204 }
2205 
2206 
2207 void Foam::conformalVoronoiMesh::storeSurfaceConformation()
2208 {
2209  Info<< nl << "Storing surface conformation" << endl;
2210 
2211  surfaceConformationVertices_.clear();
2212 
2213  // Use a temporary dynamic list to speed up insertion.
2214  DynamicList<Vb> tempSurfaceVertices(number_of_vertices()/10);
2215 
2216  for
2217  (
2218  Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
2219  vit != finite_vertices_end();
2220  vit++
2221  )
2222  {
2223  // Store points that are not referred, part of a pair, but not feature
2224  // points
2225  if
2226  (
2227  !vit->referred()
2228  && vit->boundaryPoint()
2229  && !vit->featurePoint()
2230  && !vit->constrained()
2231  )
2232  {
2233  tempSurfaceVertices.append
2234  (
2235  Vb
2236  (
2237  vit->point(),
2238  vit->index(),
2239  vit->type(),
2241  )
2242  );
2243  }
2244  }
2245 
2246  tempSurfaceVertices.shrink();
2247 
2248  surfaceConformationVertices_.transfer(tempSurfaceVertices);
2249 
2250  Info<< " Stored "
2251  << returnReduce
2252  (
2253  label(surfaceConformationVertices_.size()),
2254  sumOp<label>()
2255  )
2256  << " vertices" << nl << endl;
2257 }
2258 
2259 
2260 void Foam::conformalVoronoiMesh::reinsertSurfaceConformation()
2261 {
2262  Info<< nl << "Reinserting stored surface conformation" << endl;
2263 
2264  Map<label> oldToNewIndices =
2265  insertPointPairs(surfaceConformationVertices_, true, true);
2266 
2267  ptPairs_.reIndex(oldToNewIndices);
2268 
2269  bitSet selectedElems(surfaceConformationVertices_.size(), true);
2270 
2271  forAll(surfaceConformationVertices_, vI)
2272  {
2273  Vb& v = surfaceConformationVertices_[vI];
2274  label& vIndex = v.index();
2275 
2276  const auto iter = oldToNewIndices.cfind(vIndex);
2277 
2278  if (iter.found())
2279  {
2280  const label newIndex = *iter;
2281 
2282  if (newIndex != -1)
2283  {
2284  vIndex = newIndex;
2285  }
2286  else
2287  {
2288  selectedElems.unset(vI);
2289  }
2290  }
2291  }
2292 
2293  inplaceSubset<bitSet, List<Vb>>
2294  (
2295  selectedElems,
2296  surfaceConformationVertices_
2297  );
2298 }
2299 
2300 
2301 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::uniform
static Vector< Cmpt > uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:164
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Definition: gatherScatterList.C:215
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::vectorTools::cosPhi
T cosPhi(const Vector< T > &a, const Vector< T > &b, const T &tolerance=SMALL)
Calculate angle between a and b in radians.
Definition: vectorTools.H:107
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
backgroundMeshDecomposition.H
Foam::constant::physicoChemical::c1
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
CGAL::indexedVertex::index
Foam::label & index()
Definition: indexedVertexI.H:159
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
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::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::vectorTools
Collection of functions for testing relationships between two vectors.
Definition: vectorTools.H:49
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
Foam::UPstream::allProcs
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition: UPstream.H:508
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Vector< scalar >
Foam::labelPairHashSet
HashSet< labelPair, Foam::Hash< labelPair > > labelPairHashSet
A HashSet for a labelPair. The hashing is based on labelPair (FixedList) and is thus non-commutative.
Definition: labelPairHashes.H:65
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
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::limit
complex limit(const complex &, const complex &)
Definition: complexI.H:263
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
vectorTools.H
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
OBJstream.H
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
conformalVoronoiMesh.H
indexedCellChecks.H
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265