DistributedDelaunayMesh.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 
30 #include "meshSearch.H"
31 #include "mapDistribute.H"
33 #include "pointConversion.H"
34 #include "indexedVertexEnum.H"
35 #include "IOmanip.H"
36 #include <algorithm>
37 #include <random>
38 
39 // * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * * //
40 
41 template<class Triangulation>
44 (
45  const List<label>& toProc
46 )
47 {
48  // Determine send map
49  // ~~~~~~~~~~~~~~~~~~
50 
51  // 1. Count
52  labelList nSend(Pstream::nProcs(), Zero);
53 
54  forAll(toProc, i)
55  {
56  label proci = toProc[i];
57 
58  nSend[proci]++;
59  }
60 
61 
62  // 2. Size sendMap
63  labelListList sendMap(Pstream::nProcs());
64 
65  forAll(nSend, proci)
66  {
67  sendMap[proci].setSize(nSend[proci]);
68 
69  nSend[proci] = 0;
70  }
71 
72  // 3. Fill sendMap
73  forAll(toProc, i)
74  {
75  label proci = toProc[i];
76 
77  sendMap[proci][nSend[proci]++] = i;
78  }
79 
80  // 4. Send over how many I need to receive
81  labelList recvSizes;
82  Pstream::exchangeSizes(sendMap, recvSizes);
83 
84 
85  // Determine receive map
86  // ~~~~~~~~~~~~~~~~~~~~~
87 
88  labelListList constructMap(Pstream::nProcs());
89 
90  // Local transfers first
91  constructMap[Pstream::myProcNo()] = identity
92  (
93  sendMap[Pstream::myProcNo()].size()
94  );
95 
96  label constructSize = constructMap[Pstream::myProcNo()].size();
97 
98  forAll(constructMap, proci)
99  {
100  if (proci != Pstream::myProcNo())
101  {
102  label nRecv = recvSizes[proci];
103 
104  constructMap[proci].setSize(nRecv);
105 
106  for (label i = 0; i < nRecv; i++)
107  {
108  constructMap[proci][i] = constructSize++;
109  }
110  }
111  }
112 
114  (
115  constructSize,
116  std::move(sendMap),
117  std::move(constructMap)
118  );
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123 
124 template<class Triangulation>
126 (
127  const Time& runTime
128 )
129 :
130  DelaunayMesh<Triangulation>(runTime),
131  allBackgroundMeshBounds_()
132 {}
133 
134 
135 template<class Triangulation>
137 (
138  const Time& runTime,
139  const word& meshName
140 )
141 :
142  DelaunayMesh<Triangulation>(runTime, meshName),
143  allBackgroundMeshBounds_()
144 {}
145 
146 
147 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
148 
149 template<class Triangulation>
151 (
152  const boundBox& bb
153 )
154 {
155  allBackgroundMeshBounds_.reset(new List<boundBox>(Pstream::nProcs()));
156 
157  // Give the bounds of every processor to every other processor
158  allBackgroundMeshBounds_()[Pstream::myProcNo()] = bb;
159 
160  Pstream::gatherList(allBackgroundMeshBounds_());
161  Pstream::scatterList(allBackgroundMeshBounds_());
162 
163  return true;
164 }
165 
166 
167 template<class Triangulation>
169 (
170  const Vertex_handle& v
171 ) const
172 {
173  return isLocal(v->procIndex());
174 }
175 
176 
177 template<class Triangulation>
179 (
180  const label localProcIndex
181 ) const
182 {
183  return localProcIndex == Pstream::myProcNo();
184 }
185 
186 
187 template<class Triangulation>
189 (
190  const point& centre,
191  const scalar radiusSqr
192 ) const
193 {
194  DynamicList<label> toProc(Pstream::nProcs());
195 
196  forAll(allBackgroundMeshBounds_(), proci)
197  {
198  // Test against the bounding box of the processor
199  if
200  (
201  !isLocal(proci)
202  && allBackgroundMeshBounds_()[proci].overlaps(centre, radiusSqr)
203  )
204  {
205  toProc.append(proci);
206  }
207  }
208 
209  return toProc;
210 }
211 
212 
213 template<class Triangulation>
215 (
216  const Cell_handle& cit,
217  Map<labelList>& circumsphereOverlaps
218 ) const
219 {
220  const Foam::point& cc = cit->dual();
221 
222  const scalar crSqr = magSqr
223  (
224  cc - topoint(cit->vertex(0)->point())
225  );
226 
227  labelList circumsphereOverlap = overlapProcessors
228  (
229  cc,
230  sqr(1.01)*crSqr
231  );
232 
233  cit->cellIndex() = this->getNewCellIndex();
234 
235  if (!circumsphereOverlap.empty())
236  {
237  circumsphereOverlaps.insert(cit->cellIndex(), circumsphereOverlap);
238 
239  return true;
240  }
241 
242  return false;
243 }
244 
245 
246 template<class Triangulation>
248 (
249  Map<labelList>& circumsphereOverlaps
250 ) const
251 {
252  // Start by assuming that all the cells have no index
253  // If they do, they have already been visited so ignore them
254 
255  labelHashSet cellToCheck
256  (
257  Triangulation::number_of_finite_cells()
258  /Pstream::nProcs()
259  );
260 
261 // std::list<Cell_handle> infinite_cells;
262 // Triangulation::incident_cells
263 // (
264 // Triangulation::infinite_vertex(),
265 // std::back_inserter(infinite_cells)
266 // );
267 //
268 // for
269 // (
270 // typename std::list<Cell_handle>::iterator vcit
271 // = infinite_cells.begin();
272 // vcit != infinite_cells.end();
273 // ++vcit
274 // )
275 // {
276 // Cell_handle cit = *vcit;
277 //
278 // // Index of infinite vertex in this cell.
279 // label i = cit->index(Triangulation::infinite_vertex());
280 //
281 // Cell_handle c = cit->neighbor(i);
282 //
283 // if (c->unassigned())
284 // {
285 // c->cellIndex() = this->getNewCellIndex();
286 //
287 // if (checkProcBoundaryCell(c, circumsphereOverlaps))
288 // {
289 // cellToCheck.insert(c->cellIndex());
290 // }
291 // }
292 // }
293 //
294 //
295 // for
296 // (
297 // Finite_cells_iterator cit = Triangulation::finite_cells_begin();
298 // cit != Triangulation::finite_cells_end();
299 // ++cit
300 // )
301 // {
302 // if (cit->parallelDualVertex())
303 // {
304 // if (cit->unassigned())
305 // {
306 // if (checkProcBoundaryCell(cit, circumsphereOverlaps))
307 // {
308 // cellToCheck.insert(cit->cellIndex());
309 // }
310 // }
311 // }
312 // }
313 
314 
315  for
316  (
317  All_cells_iterator cit = Triangulation::all_cells_begin();
318  cit != Triangulation::all_cells_end();
319  ++cit
320  )
321  {
322  if (Triangulation::is_infinite(cit))
323  {
324  // Index of infinite vertex in this cell.
325  label i = cit->index(Triangulation::infinite_vertex());
326 
327  Cell_handle c = cit->neighbor(i);
328 
329  if (c->unassigned())
330  {
331  c->cellIndex() = this->getNewCellIndex();
332 
333  if (checkProcBoundaryCell(c, circumsphereOverlaps))
334  {
335  cellToCheck.insert(c->cellIndex());
336  }
337  }
338  }
339  else if (cit->parallelDualVertex())
340  {
341  if (cit->unassigned())
342  {
343  if (checkProcBoundaryCell(cit, circumsphereOverlaps))
344  {
345  cellToCheck.insert(cit->cellIndex());
346  }
347  }
348  }
349  }
350 
351  for
352  (
353  Finite_cells_iterator cit = Triangulation::finite_cells_begin();
354  cit != Triangulation::finite_cells_end();
355  ++cit
356  )
357  {
358  if (cellToCheck.found(cit->cellIndex()))
359  {
360  // Get the neighbours and check them
361  for (label adjCelli = 0; adjCelli < 4; ++adjCelli)
362  {
363  Cell_handle citNeighbor = cit->neighbor(adjCelli);
364 
365  // Ignore if has far point or previously visited
366  if
367  (
368  !citNeighbor->unassigned()
369  || !citNeighbor->internalOrBoundaryDualVertex()
370  || Triangulation::is_infinite(citNeighbor)
371  )
372  {
373  continue;
374  }
375 
376  if
377  (
378  checkProcBoundaryCell
379  (
380  citNeighbor,
381  circumsphereOverlaps
382  )
383  )
384  {
385  cellToCheck.insert(citNeighbor->cellIndex());
386  }
387  }
388 
389  cellToCheck.unset(cit->cellIndex());
390  }
391  }
392 }
393 
394 
395 template<class Triangulation>
397 (
398  const Map<labelList>& circumsphereOverlaps,
399  PtrList<labelPairHashSet>& referralVertices,
400  DynamicList<label>& targetProcessor,
401  DynamicList<Vb>& parallelInfluenceVertices
402 )
403 {
404  // Relying on the order of iteration of cells being the same as before
405  for
406  (
407  Finite_cells_iterator cit = Triangulation::finite_cells_begin();
408  cit != Triangulation::finite_cells_end();
409  ++cit
410  )
411  {
412  if (Triangulation::is_infinite(cit))
413  {
414  continue;
415  }
416 
417  const auto iter = circumsphereOverlaps.cfind(cit->cellIndex());
418 
419  // Pre-tested circumsphere potential influence
420  if (iter.found())
421  {
422  const labelList& citOverlaps = iter();
423 
424  for (const label proci : citOverlaps)
425  {
426  for (int i = 0; i < 4; i++)
427  {
428  Vertex_handle v = cit->vertex(i);
429 
430  if (v->farPoint())
431  {
432  continue;
433  }
434 
435  label vProcIndex = v->procIndex();
436  label vIndex = v->index();
437 
438  const labelPair procIndexPair(vProcIndex, vIndex);
439 
440  // Using the hashSet to ensure that each vertex is only
441  // referred once to each processor.
442  // Do not refer a vertex to its own processor.
443  if (vProcIndex != proci)
444  {
445  if (referralVertices[proci].insert(procIndexPair))
446  {
447  targetProcessor.append(proci);
448 
449  parallelInfluenceVertices.append
450  (
451  Vb
452  (
453  v->point(),
454  v->index(),
455  v->type(),
456  v->procIndex()
457  )
458  );
459 
460  parallelInfluenceVertices.last().targetCellSize() =
461  v->targetCellSize();
462  parallelInfluenceVertices.last().alignment() =
463  v->alignment();
464  }
465  }
466  }
467  }
468  }
469  }
470 }
471 
472 
473 template<class Triangulation>
475 (
476  const DynamicList<label>& targetProcessor,
477  DynamicList<Vb>& parallelVertices,
478  PtrList<labelPairHashSet>& referralVertices,
479  labelPairHashSet& receivedVertices
480 )
481 {
482  DynamicList<Vb> referredVertices(targetProcessor.size());
483 
484  const label preDistributionSize = parallelVertices.size();
485 
486  autoPtr<mapDistribute> pointMapPtr = buildMap(targetProcessor);
487  mapDistribute& pointMap = *pointMapPtr;
488 
489  // Make a copy of the original list.
490  DynamicList<Vb> originalParallelVertices(parallelVertices);
491 
492  pointMap.distribute(parallelVertices);
493 
494  for (const int proci : Pstream::allProcs())
495  {
496  const labelList& constructMap = pointMap.constructMap()[proci];
497 
498  if (constructMap.size())
499  {
500  forAll(constructMap, i)
501  {
502  const Vb& v = parallelVertices[constructMap[i]];
503 
504  if
505  (
506  v.procIndex() != Pstream::myProcNo()
507  && !receivedVertices.found(labelPair(v.procIndex(), v.index()))
508  )
509  {
510  referredVertices.append(v);
511 
512  receivedVertices.insert
513  (
514  labelPair(v.procIndex(), v.index())
515  );
516  }
517  }
518  }
519  }
520 
521  label preInsertionSize = Triangulation::number_of_vertices();
522 
523  labelPairHashSet pointsNotInserted = rangeInsertReferredWithInfo
524  (
525  referredVertices.begin(),
526  referredVertices.end(),
527  true
528  );
529 
530  if (!pointsNotInserted.empty())
531  {
532  forAllConstIters(pointsNotInserted, iter)
533  {
534  if (receivedVertices.found(iter.key()))
535  {
536  receivedVertices.erase(iter.key());
537  }
538  }
539  }
540 
541  boolList pointInserted(parallelVertices.size(), true);
542 
543  forAll(parallelVertices, vI)
544  {
545  const labelPair procIndexI
546  (
547  parallelVertices[vI].procIndex(),
548  parallelVertices[vI].index()
549  );
550 
551  if (pointsNotInserted.found(procIndexI))
552  {
553  pointInserted[vI] = false;
554  }
555  }
556 
557  pointMap.reverseDistribute(preDistributionSize, pointInserted);
558 
559  forAll(originalParallelVertices, vI)
560  {
561  const label procIndex = targetProcessor[vI];
562 
563  if (!pointInserted[vI])
564  {
565  if (referralVertices[procIndex].size())
566  {
567  if
568  (
569  !referralVertices[procIndex].unset
570  (
571  labelPair
572  (
573  originalParallelVertices[vI].procIndex(),
574  originalParallelVertices[vI].index()
575  )
576  )
577  )
578  {
579  Pout<< "*** not found "
580  << originalParallelVertices[vI].procIndex()
581  << " " << originalParallelVertices[vI].index() << endl;
582  }
583  }
584  }
585  }
586 
587  label postInsertionSize = Triangulation::number_of_vertices();
588 
589  reduce(preInsertionSize, sumOp<label>());
590  reduce(postInsertionSize, sumOp<label>());
591 
592  label nTotalToInsert = referredVertices.size();
593 
594  reduce(nTotalToInsert, sumOp<label>());
595 
596  if (preInsertionSize + nTotalToInsert != postInsertionSize)
597  {
598  label nNotInserted =
599  returnReduce(pointsNotInserted.size(), sumOp<label>());
600 
601  Info<< " Inserted = "
602  << setw(name(label(Triangulation::number_of_finite_cells())).size())
603  << nTotalToInsert - nNotInserted
604  << " / " << nTotalToInsert << endl;
605 
606  nTotalToInsert -= nNotInserted;
607  }
608  else
609  {
610  Info<< " Inserted = " << nTotalToInsert << endl;
611  }
612 
613  return nTotalToInsert;
614 }
615 
616 
617 template<class Triangulation>
619 (
620  const boundBox& bb,
621  PtrList<labelPairHashSet>& referralVertices,
622  labelPairHashSet& receivedVertices,
623  bool iterateReferral
624 )
625 {
626  if (!Pstream::parRun())
627  {
628  return;
629  }
630 
631  if (!allBackgroundMeshBounds_)
632  {
633  distributeBoundBoxes(bb);
634  }
635 
636  label nVerts = Triangulation::number_of_vertices();
637  label nCells = Triangulation::number_of_finite_cells();
638 
639  DynamicList<Vb> parallelInfluenceVertices(0.1*nVerts);
640  DynamicList<label> targetProcessor(0.1*nVerts);
641 
642  // Some of these values will not be used, i.e. for non-real cells
643  DynamicList<Foam::point> circumcentre(0.1*nVerts);
644  DynamicList<scalar> circumradiusSqr(0.1*nVerts);
645 
646  Map<labelList> circumsphereOverlaps(nCells);
647 
648  findProcessorBoundaryCells(circumsphereOverlaps);
649 
650  Info<< " Influences = "
651  << setw(name(nCells).size())
652  << returnReduce(circumsphereOverlaps.size(), sumOp<label>()) << " / "
653  << returnReduce(nCells, sumOp<label>());
654 
655  markVerticesToRefer
656  (
657  circumsphereOverlaps,
658  referralVertices,
659  targetProcessor,
660  parallelInfluenceVertices
661  );
662 
663  referVertices
664  (
665  targetProcessor,
666  parallelInfluenceVertices,
667  referralVertices,
668  receivedVertices
669  );
670 
671  if (iterateReferral)
672  {
673  label oldNReferred = 0;
674  label nIterations = 1;
675 
676  Info<< incrIndent << indent
677  << "Iteratively referring referred vertices..."
678  << endl;
679  do
680  {
681  Info<< indent << "Iteration " << nIterations++ << ":";
682 
683  circumsphereOverlaps.clear();
684  targetProcessor.clear();
685  parallelInfluenceVertices.clear();
686 
687  findProcessorBoundaryCells(circumsphereOverlaps);
688 
689  nCells = Triangulation::number_of_finite_cells();
690 
691  Info<< " Influences = "
692  << setw(name(nCells).size())
693  << returnReduce(circumsphereOverlaps.size(), sumOp<label>())
694  << " / "
695  << returnReduce(nCells, sumOp<label>());
696 
697  markVerticesToRefer
698  (
699  circumsphereOverlaps,
700  referralVertices,
701  targetProcessor,
702  parallelInfluenceVertices
703  );
704 
705  label nReferred = referVertices
706  (
707  targetProcessor,
708  parallelInfluenceVertices,
709  referralVertices,
710  receivedVertices
711  );
712 
713  if (nReferred == 0 || nReferred == oldNReferred)
714  {
715  break;
716  }
717 
718  oldNReferred = nReferred;
719 
720  } while (true);
721 
722  Info<< decrIndent;
723  }
724 }
725 
726 
727 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
728 
729 
730 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
731 
732 template<class Triangulation>
733 Foam::scalar
735 {
736  label nRealVertices = 0;
737 
738  for
739  (
740  Finite_vertices_iterator vit = Triangulation::finite_vertices_begin();
741  vit != Triangulation::finite_vertices_end();
742  ++vit
743  )
744  {
745  // Only store real vertices that are not feature vertices
746  if (vit->real() && !vit->featurePoint())
747  {
748  nRealVertices++;
749  }
750  }
751 
752  scalar globalNRealVertices = returnReduce
753  (
754  nRealVertices,
755  sumOp<label>()
756  );
757 
758  scalar unbalance = returnReduce
759  (
760  mag(1.0 - nRealVertices/(globalNRealVertices/Pstream::nProcs())),
761  maxOp<scalar>()
762  );
763 
764  Info<< " Processor unbalance " << unbalance << endl;
765 
766  return unbalance;
767 }
768 
769 
770 template<class Triangulation>
772 (
773  const boundBox& bb
774 )
775 {
777 
778  if (!Pstream::parRun())
779  {
780  return false;
781  }
782 
783  distributeBoundBoxes(bb);
784 
785  return true;
786 }
787 
788 
789 template<class Triangulation>
792 (
793  const backgroundMeshDecomposition& decomposition,
794  List<Foam::point>& points
795 )
796 {
797  if (!Pstream::parRun())
798  {
799  return nullptr;
800  }
801 
802  distributeBoundBoxes(decomposition.procBounds());
803 
804  return decomposition.distributePoints(points);
805 }
806 
807 
808 template<class Triangulation>
810 {
811  if (!Pstream::parRun())
812  {
813  return;
814  }
815 
816  if (!allBackgroundMeshBounds_)
817  {
818  distributeBoundBoxes(bb);
819  }
820 
821  const label nApproxReferred =
822  Triangulation::number_of_vertices()
823  /Pstream::nProcs();
824 
825  PtrList<labelPairHashSet> referralVertices(Pstream::nProcs());
826  forAll(referralVertices, proci)
827  {
828  if (!isLocal(proci))
829  {
830  referralVertices.set(proci, new labelPairHashSet(nApproxReferred));
831  }
832  }
833 
834  labelPairHashSet receivedVertices(nApproxReferred);
835 
836  sync
837  (
838  bb,
839  referralVertices,
840  receivedVertices,
841  true
842  );
843 }
844 
845 
846 template<class Triangulation>
847 template<class PointIterator>
850 (
851  PointIterator begin,
852  PointIterator end,
853  bool printErrors
854 )
855 {
856  const boundBox& bb = allBackgroundMeshBounds_()[Pstream::myProcNo()];
857 
858  typedef DynamicList
859  <
860  std::pair<scalar, label>
861  > vectorPairPointIndex;
862 
863  vectorPairPointIndex pointsBbDistSqr;
864 
865  label count = 0;
866  for (PointIterator it = begin; it != end; ++it)
867  {
868  const Foam::point samplePoint(topoint(it->point()));
869 
870  scalar distFromBbSqr = 0;
871 
872  if (!bb.contains(samplePoint))
873  {
874  const Foam::point nearestPoint = bb.nearest(samplePoint);
875 
876  distFromBbSqr = magSqr(nearestPoint - samplePoint);
877  }
878 
879  pointsBbDistSqr.append
880  (
881  std::make_pair(distFromBbSqr, count++)
882  );
883  }
884 
886  (
887  pointsBbDistSqr.begin(),
888  pointsBbDistSqr.end(),
889  std::default_random_engine()
890  );
891 
892  // Sort in ascending order by the distance of the point from the centre
893  // of the processor bounding box
894  sort(pointsBbDistSqr.begin(), pointsBbDistSqr.end());
895 
896  typename Triangulation::Vertex_handle hint;
897 
898  typename Triangulation::Locate_type lt;
899  int li, lj;
900 
901  label nNotInserted = 0;
902 
903  labelPairHashSet uninserted
904  (
905  Triangulation::number_of_vertices()
906  /Pstream::nProcs()
907  );
908 
909  for
910  (
911  typename vectorPairPointIndex::const_iterator p =
912  pointsBbDistSqr.begin();
913  p != pointsBbDistSqr.end();
914  ++p
915  )
916  {
917  const size_t checkInsertion = Triangulation::number_of_vertices();
918 
919  const Vb& vert = *(begin + p->second);
920  const Point& pointToInsert = vert.point();
921 
922  // Locate the point
923  Cell_handle c = Triangulation::locate(pointToInsert, lt, li, lj, hint);
924 
925  bool inserted = false;
926 
927  if (lt == Triangulation::VERTEX)
928  {
929  if (printErrors)
930  {
931  Vertex_handle nearV =
932  Triangulation::nearest_vertex(pointToInsert);
933 
934  Pout<< "Failed insertion, point already exists" << nl
935  << "Failed insertion : " << vert.info()
936  << " nearest : " << nearV->info();
937  }
938  }
939  else if (lt == Triangulation::OUTSIDE_AFFINE_HULL)
940  {
942  << "Point is outside affine hull! pt = " << pointToInsert
943  << endl;
944  }
945  else if (lt == Triangulation::OUTSIDE_CONVEX_HULL)
946  {
947  // TODO: Can this be optimised?
948  //
949  // Only want to insert if a connection is formed between
950  // pointToInsert and an internal or internal boundary point.
951  hint = Triangulation::insert(pointToInsert, c);
952  inserted = true;
953  }
954  else
955  {
956  // Get the cells that conflict with p in a vector V,
957  // and a facet on the boundary of this hole in f.
958  std::vector<Cell_handle> V;
959  typename Triangulation::Facet f;
960 
961  Triangulation::find_conflicts
962  (
963  pointToInsert,
964  c,
965  CGAL::Oneset_iterator<typename Triangulation::Facet>(f),
966  std::back_inserter(V)
967  );
968 
969  for (size_t i = 0; i < V.size(); ++i)
970  {
971  Cell_handle conflictingCell = V[i];
972 
973  if
974  (
975  Triangulation::dimension() < 3 // 2D triangulation
976  ||
977  (
978  !Triangulation::is_infinite(conflictingCell)
979  && (
980  conflictingCell->real()
981  || conflictingCell->hasFarPoint()
982  )
983  )
984  )
985  {
986  hint = Triangulation::insert_in_hole
987  (
988  pointToInsert,
989  V.begin(),
990  V.end(),
991  f.first,
992  f.second
993  );
994 
995  inserted = true;
996 
997  break;
998  }
999  }
1000  }
1001 
1002  if (inserted)
1003  {
1004  if (checkInsertion != Triangulation::number_of_vertices() - 1)
1005  {
1006  if (printErrors)
1007  {
1008  Vertex_handle nearV =
1009  Triangulation::nearest_vertex(pointToInsert);
1010 
1011  Pout<< "Failed insertion : " << vert.info()
1012  << " nearest : " << nearV->info();
1013  }
1014  }
1015  else
1016  {
1017  hint->index() = vert.index();
1018  hint->type() = vert.type();
1019  hint->procIndex() = vert.procIndex();
1020  hint->targetCellSize() = vert.targetCellSize();
1021  hint->alignment() = vert.alignment();
1022  }
1023  }
1024  else
1025  {
1026  uninserted.insert(labelPair(vert.procIndex(), vert.index()));
1027  nNotInserted++;
1028  }
1029  }
1030 
1031  return uninserted;
1032 }
1033 
1034 
1035 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
insert
srcOptions insert("case", fileName(rootDirSource/caseDirSource))
runTime
engineTime & runTime
Definition: createEngineTime.H:13
CGAL::indexedVertex::alignment
Foam::tensor & alignment()
Definition: indexedVertexI.H:189
p
volScalarField & p
Definition: createFieldRefs.H:8
pointConversion.H
stdFoam::begin
constexpr auto begin(C &c) -> decltype(c.begin())
Return iterator to the beginning of the container c.
Definition: stdFoam.H:97
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
CGAL::indexedVertex::type
vertexType & type()
Definition: indexedVertexI.H:174
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::boolList
List< bool > boolList
A List of bools.
Definition: List.H:65
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::BitOps::unset
void unset(List< bool > &bools, const labelRange &range)
Unset the specified range 'on' in a boolList.
Definition: BitOps.C:96
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::HashSet
A HashTable with keys but without contents that is similar to std::unordered_set.
Definition: HashSet.H:77
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:346
Foam::DistributedDelaunayMesh::rangeInsertReferredWithInfo
labelPairHashSet rangeInsertReferredWithInfo(PointIterator begin, PointIterator end, bool printErrors=true)
Inserts points into the triangulation if the point is within.
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::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::labelPair
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:54
Foam::DistributedDelaunayMesh::buildMap
static autoPtr< mapDistribute > buildMap(const List< label > &toProc)
Build a mapDistribute for the supplied destination processor data.
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
indexedVertexEnum.H
Foam::DistributedDelaunayMesh::sync
void sync(const boundBox &bb)
Refer vertices so that the processor interfaces are consistent.
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
CGAL::indexedVertex::index
Foam::label & index()
Definition: indexedVertexI.H:159
Foam::sort
void sort(UList< T > &a)
Definition: UList.C:261
IOmanip.H
Istream and Ostream manipulators taking arguments.
DistributedDelaunayMesh.H
Foam::DistributedDelaunayMesh
Definition: DistributedDelaunayMesh.H:58
Foam::DistributedDelaunayMesh::calculateLoadUnbalance
scalar calculateLoadUnbalance() const
CGAL::indexedVertex::info
Foam::InfoProxy< indexedVertex< Gt, Vb > > info() const
Info proxy, to print information to a stream.
Definition: indexedVertex.H:307
CGAL::indexedVertex::targetCellSize
Foam::scalar & targetCellSize()
Definition: indexedVertexI.H:203
reduce
reduce(hasMovingMesh, orOp< bool >())
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:353
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:339
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::autoPtr< Foam::mapDistribute >
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
CGAL::indexedVertex::procIndex
int procIndex() const
Definition: indexedVertexI.H:251
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
meshSearch.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::DistributedDelaunayMesh::distribute
bool distribute(const boundBox &bb)
mapDistribute.H
f
labelList f(nPoints)
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::Vector
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:62
Foam::shuffle
void shuffle(UList< T > &a)
Definition: UList.C:289
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::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)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Point
CGAL::Point_3< K > Point
Definition: CGALIndexedPolyhedron.H:53
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
zeroGradientFvPatchFields.H
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::DelaunayMesh::reset
void reset()
Clear the entire triangulation.