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