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