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