Loading...
Searching...
No Matches
conformalVoronoiMeshCalcDualMesh.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) 2018-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 "motionSmoother.H"
32#include "polyMeshGeometry.H"
33#include "indexedCellChecks.H"
34#include "OBJstream.H"
35#include "indexedCellOps.H"
36#include "ListOps.H"
37#include "DelaunayMeshTools.H"
38
39// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
40
41void Foam::conformalVoronoiMesh::calcDualMesh
42(
43 pointField& points,
44 labelList& boundaryPts,
45 faceList& faces,
46 labelList& owner,
47 labelList& neighbour,
48 wordList& patchNames,
49 PtrList<dictionary>& patchDicts,
50 pointField& cellCentres,
51 labelList& cellToDelaunayVertex,
52 labelListList& patchToDelaunayVertex,
53 bitSet& boundaryFacesToRemove
54)
55{
56 timeCheck("Start calcDualMesh");
57
58 setVertexSizeAndAlignment();
59
60 timeCheck("After setVertexSizeAndAlignment");
61
62 indexDualVertices(points, boundaryPts);
63
64 {
65 Info<< nl << "Merging identical points" << endl;
66
67 // There is no guarantee that a merge of close points is no-risk
68 mergeIdenticalDualVertices(points, boundaryPts);
69 }
70
71 // Final dual face and owner neighbour construction
72
73 timeCheck("Before createFacesOwnerNeighbourAndPatches");
74
75 createFacesOwnerNeighbourAndPatches
76 (
77 points,
78 faces,
79 owner,
80 neighbour,
83 patchToDelaunayVertex, // from patch face to Delaunay vertex (slavePp)
84 boundaryFacesToRemove,
85 false
86 );
87
88 // deferredCollapseFaceSet(owner, neighbour, deferredCollapseFaces);
89
90 cellCentres = DelaunayMeshTools::allPoints(*this);
91
92 cellToDelaunayVertex = removeUnusedCells(owner, neighbour);
93
94 cellCentres = pointField(cellCentres, cellToDelaunayVertex);
95
96 removeUnusedPoints(faces, points, boundaryPts);
97
98 timeCheck("End of calcDualMesh");
99}
100
101
102void Foam::conformalVoronoiMesh::calcTetMesh
103(
105 labelList& pointToDelaunayVertex,
106 faceList& faces,
107 labelList& owner,
108 labelList& neighbour,
111)
112{
113 labelList vertexMap(number_of_vertices());
114
115 label vertI = 0;
116
117 points.setSize(number_of_vertices());
118 pointToDelaunayVertex.setSize(number_of_vertices());
119
120 for
121 (
122 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
123 vit != finite_vertices_end();
124 ++vit
125 )
126 {
127 if (vit->internalPoint() || vit->boundaryPoint())
128 {
129 vertexMap[vit->index()] = vertI;
130 points[vertI] = topoint(vit->point());
131 pointToDelaunayVertex[vertI] = vit->index();
132 vertI++;
133 }
134 }
135
136 points.setSize(vertI);
137 pointToDelaunayVertex.setSize(vertI);
138
139 label celli = 0;
140
141 for
142 (
143 Delaunay::Finite_cells_iterator cit = finite_cells_begin();
144 cit != finite_cells_end();
145 ++cit
146 )
147 {
148 if (cit->internalOrBoundaryDualVertex())
149 {
150 cit->cellIndex() = celli++;
151 }
152 else
153 {
154 cit->cellIndex() = Cb::ctFar;
155 }
156 }
157
158 patchNames = geometryToConformTo_.patchNames();
159
161
162 patchNames[patchNames.size() - 1] = "foamyHexMesh_defaultPatch";
163
164 label nPatches = patchNames.size();
165
166 List<DynamicList<face>> patchFaces(nPatches);
167 List<DynamicList<label>> patchOwners(nPatches);
168
169 faces.setSize(number_of_finite_facets());
170
171 owner.setSize(number_of_finite_facets());
172
173 neighbour.setSize(number_of_finite_facets());
174
175 label facei = 0;
176
177 labelList verticesOnTriFace(3, label(-1));
178
179 face newFace(verticesOnTriFace);
180
181 for
182 (
183 Delaunay::Finite_facets_iterator fit = finite_facets_begin();
184 fit != finite_facets_end();
185 ++fit
186 )
187 {
188 const Cell_handle c1(fit->first);
189 const label oppositeVertex = fit->second;
190 const Cell_handle c2(c1->neighbor(oppositeVertex));
191
192 if (c1->hasFarPoint() && c2->hasFarPoint())
193 {
194 // Both tets are outside, skip
195 continue;
196 }
197
198 label c1I = c1->cellIndex();
199 label c2I = c2->cellIndex();
200
201 label ownerCell = -1;
202 label neighbourCell = -1;
203
204 for (label i = 0; i < 3; i++)
205 {
206 verticesOnTriFace[i] = vertexMap
207 [
208 c1->vertex(vertex_triple_index(oppositeVertex, i))->index()
209 ];
210 }
211
212 newFace = face(verticesOnTriFace);
213
214 if (c1->hasFarPoint() || c2->hasFarPoint())
215 {
216 // Boundary face...
217 if (c1->hasFarPoint())
218 {
219 //... with c1 outside
220 ownerCell = c2I;
221 }
222 else
223 {
224 // ... with c2 outside
225 ownerCell = c1I;
226
227 reverse(newFace);
228 }
229
230 label patchIndex = geometryToConformTo_.findPatch
231 (
232 newFace.centre(points)
233 );
234
235 if (patchIndex == -1)
236 {
237 patchIndex = patchNames.size() - 1;
238
240 << newFace.centre(points) << nl
241 << "did not find a surface patch. Adding to "
242 << patchNames[patchIndex]
243 << endl;
244 }
245
246 patchFaces[patchIndex].append(newFace);
247 patchOwners[patchIndex].append(ownerCell);
248 }
249 else
250 {
251 // Internal face...
252 if (c1I < c2I)
253 {
254 // ...with c1 as the ownerCell
255 ownerCell = c1I;
256 neighbourCell = c2I;
257
258 reverse(newFace);
259 }
260 else
261 {
262 // ...with c2 as the ownerCell
263 ownerCell = c2I;
264 neighbourCell = c1I;
265 }
266
267 faces[facei] = newFace;
268 owner[facei] = ownerCell;
269 neighbour[facei] = neighbourCell;
270 facei++;
271 }
272 }
273
274 label nInternalFaces = facei;
275
276 faces.setSize(nInternalFaces);
277 owner.setSize(nInternalFaces);
278 neighbour.setSize(nInternalFaces);
279
280 sortFaces(faces, owner, neighbour);
281
282// bitSet boundaryFacesToRemove;
283// List<DynamicList<bool>> indirectPatchFace;
284//
285// addPatches
286// (
287// nInternalFaces,
288// faces,
289// owner,
290// patchDicts,
291// boundaryFacesToRemove,
292// patchFaces,
293// patchOwners,
294// indirectPatchFace
295// );
296}
297
298
299void Foam::conformalVoronoiMesh::mergeIdenticalDualVertices
300(
301 const pointField& pts,
302 labelList& boundaryPts
303)
304{
305 // Assess close points to be merged
306
307 label nPtsMerged = 0;
308 label nPtsMergedSum = 0;
309
310 do
311 {
312 Map<label> dualPtIndexMap;
313
314 nPtsMerged = mergeIdenticalDualVertices
315 (
316 pts,
317 dualPtIndexMap
318 );
319
320 reindexDualVertices(dualPtIndexMap, boundaryPts);
321
322 reduce(nPtsMerged, sumOp<label>());
323
324 nPtsMergedSum += nPtsMerged;
325
326 } while (nPtsMerged > 0);
327
328 if (nPtsMergedSum > 0)
329 {
330 Info<< " Merged " << nPtsMergedSum << " points " << endl;
331 }
332}
333
334
335Foam::label Foam::conformalVoronoiMesh::mergeIdenticalDualVertices
336(
337 const pointField& pts,
338 Map<label>& dualPtIndexMap
339) const
340{
341 label nPtsMerged = 0;
342
343 for
344 (
345 Delaunay::Finite_facets_iterator fit = finite_facets_begin();
346 fit != finite_facets_end();
347 ++fit
348 )
349 {
350 const Cell_handle c1(fit->first);
351 const label oppositeVertex = fit->second;
352 const Cell_handle c2(c1->neighbor(oppositeVertex));
353
354 if (is_infinite(c1) || is_infinite(c2))
355 {
356 continue;
357 }
358
359 label& c1I = c1->cellIndex();
360 label& c2I = c2->cellIndex();
361
362 if ((c1I != c2I) && !c1->hasFarPoint() && !c2->hasFarPoint())
363 {
364 const Foam::point& p1 = pts[c1I];
365 const Foam::point& p2 = pts[c2I];
366
367 if (p1 == p2)
368 {
369// if (c1->parallelDualVertex() || c2->parallelDualVertex())
370// {
371// if (c1->vertexLowestProc() < c2->vertexLowestProc())
372// {
373// dualPtIndexMap.insert(c1I, c1I);
374// dualPtIndexMap.insert(c2I, c1I);
375// }
376// else
377// {
378// dualPtIndexMap.insert(c1I, c2I);
379// dualPtIndexMap.insert(c2I, c2I);
380// }
381// }
382 if (c1I < c2I)
383 {
384 dualPtIndexMap.insert(c1I, c1I);
385 dualPtIndexMap.insert(c2I, c1I);
386 }
387 else
388 {
389 dualPtIndexMap.insert(c1I, c2I);
390 dualPtIndexMap.insert(c2I, c2I);
391 }
392
393 nPtsMerged++;
394 }
395 }
396 }
397
398 if (debug)
399 {
400 Info<< "mergeIdenticalDualVertices:" << endl
401 << " zero-length edges : "
402 << returnReduce(nPtsMerged, sumOp<label>()) << endl
403 << endl;
404 }
405
406 return nPtsMerged;
407}
408
409
410//void Foam::conformalVoronoiMesh::smoothSurface
411//(
412// pointField& pts,
413// const labelList& boundaryPts
414//)
415//{
416// label nCollapsedFaces = 0;
417//
418// label iterI = 0;
419//
420// do
421// {
422// Map<label> dualPtIndexMap;
423//
424// nCollapsedFaces = smoothSurfaceDualFaces
425// (
426// pts,
427// boundaryPts,
428// dualPtIndexMap
429// );
430//
431// reduce(nCollapsedFaces, sumOp<label>());
432//
433// reindexDualVertices(dualPtIndexMap);
434//
435// mergeIdenticalDualVertices(pts, boundaryPts);
436//
437// if (nCollapsedFaces > 0)
438// {
439// Info<< " Collapsed " << nCollapsedFaces << " boundary faces"
440// << endl;
441// }
442//
443// if (++iterI > foamyHexMeshControls().maxCollapseIterations())
444// {
445// Info<< " maxCollapseIterations reached, stopping collapse"
446// << endl;
447//
448// break;
449// }
450//
451// } while (nCollapsedFaces > 0);
452//
453// // Force all points of boundary faces to be on the surface
504//}
505//
506//
507//Foam::label Foam::conformalVoronoiMesh::smoothSurfaceDualFaces
508//(
509// pointField& pts,
510// const labelList& boundaryPts,
511// Map<label>& dualPtIndexMap
512//) const
513//{
514// label nCollapsedFaces = 0;
515//
516// const scalar cosPerpendicularToleranceAngle = cos
517// (
518// degToRad(foamyHexMeshControls().surfaceStepFaceAngle())
519// );
520//
521// for
522// (
523// Delaunay::Finite_edges_iterator eit = finite_edges_begin();
524// eit != finite_edges_end();
525// ++eit
526// )
527// {
528// Cell_circulator ccStart = incident_cells(*eit);
529// Cell_circulator cc = ccStart;
530//
531// bool skipFace = false;
532//
533// do
534// {
535// if (dualPtIndexMap.found(cc->cellIndex()))
536// {
537// // One of the points of this face has already been
538// // collapsed this sweep, leave for next sweep
539//
540// skipFace = true;
541//
542// break;
543// }
544//
545// } while (++cc != ccStart);
546//
547// if (skipFace)
548// {
549// continue;
550// }
551//
552// if (isBoundaryDualFace(eit))
553// {
554// face dualFace = buildDualFace(eit);
555//
556// if (dualFace.size() < 3)
557// {
558// // This face has been collapsed already
559// continue;
560// }
561//
562// label maxFC = maxFilterCount(eit);
563//
564// if (maxFC > foamyHexMeshControls().filterCountSkipThreshold())
565// {
566// // A vertex on this face has been limited too many
567// // times, skip
568// continue;
569// }
570//
571// Cell_handle c = eit->first;
572// Vertex_handle vA = c->vertex(eit->second);
573// Vertex_handle vB = c->vertex(eit->third);
574//
575// if
576// (
577// vA->internalBoundaryPoint() && vA->surfacePoint()
578// && vB->externalBoundaryPoint() && vB->surfacePoint()
579// )
580// {
581// if (vA->index() == vB->index() - 1)
582// {
583// continue;
584// }
585// }
586// else if
587// (
588// vA->externalBoundaryPoint() && vA->surfacePoint()
589// && vB->internalBoundaryPoint() && vB->surfacePoint()
590// )
591// {
592// if (vA->index() == vB->index() + 1)
593// {
594// continue;
595// }
596// }
641//
642//
643// if ((faceNormal & surfaceNormal) < cosPerpendicularToleranceAngle)
644// {
645// scalar targetFaceSize = averageAnyCellSize(vA, vB);
646//
647// // Selecting faces to collapse based on angle to
648// // surface, so set collapseSizeLimitCoeff to GREAT to
649// // allow collapse of all faces
650//
651// faceCollapseMode mode = collapseFace
652// (
653// dualFace,
654// pts,
655// boundaryPts,
656// dualPtIndexMap,
657// targetFaceSize,
658// GREAT,
659// maxFC
660// );
661//
662// if (mode == fcmPoint || mode == fcmEdge)
663// {
664// nCollapsedFaces++;
665// }
666// }
667// }
668// }
669//
670// return nCollapsedFaces;
671//}
672
673
674void Foam::conformalVoronoiMesh::deferredCollapseFaceSet
675(
676 labelList& owner,
677 labelList& neighbour,
678 const labelPairHashSet& deferredCollapseFaces
679) const
680{
682
683 forAll(neighbour, nI)
684 {
685 if (deferredCollapseFaces.found(Pair<label>(owner[nI], neighbour[nI])))
686 {
687 faceLabels.append(nI);
688 }
689 }
690
691 Pout<< "facesToCollapse" << nl << faceLabels << endl;
692}
693
694
695Foam::autoPtr<Foam::polyMesh>
696Foam::conformalVoronoiMesh::createPolyMeshFromPoints
697(
698 const pointField& pts
699) const
700{
701 faceList faces;
702 labelList owner;
703 labelList neighbour;
706 pointField cellCentres;
707 labelListList patchToDelaunayVertex;
708 bitSet boundaryFacesToRemove;
709
710 timeCheck("Start of checkPolyMeshQuality");
711
712 Info<< nl << "Creating polyMesh to assess quality" << endl;
713
714 createFacesOwnerNeighbourAndPatches
715 (
716 pts,
717 faces,
718 owner,
719 neighbour,
722 patchToDelaunayVertex,
723 boundaryFacesToRemove,
724 false
725 );
726
727 cellCentres = DelaunayMeshTools::allPoints(*this);
728
729 labelList cellToDelaunayVertex(removeUnusedCells(owner, neighbour));
730 cellCentres = pointField(cellCentres, cellToDelaunayVertex);
731
733 (
735 (
736 "foamyHexMesh_temporary",
737 runTime_.timeName(),
738 runTime_,
741 ),
742 pointField(pts), // Copy of points
743 std::move(faces),
744 std::move(owner),
745 std::move(neighbour)
746 );
747
748 polyMesh& pMesh = meshPtr();
749
751
752 label nValidPatches = 0;
753
755 {
756 label nPatchFaces = patchDicts[p].get<label>("nFaces");
757
758 if
759 (
760 patchDicts.set(p)
761 && patchDicts[p].get<word>("type") == processorPolyPatch::typeName
762 )
763 {
764 // Do not create empty processor patches
765 if (nPatchFaces)
766 {
767 patchDicts[p].set("transform", "coincidentFullMatch");
768
769 patches[nValidPatches] = new processorPolyPatch
770 (
771 patchNames[p],
772 patchDicts[p],
773 nValidPatches,
774 pMesh.boundaryMesh(),
775 processorPolyPatch::typeName
776 );
777
778 nValidPatches++;
779 }
780 }
781 else
782 {
783 // Check that the patch is not empty on every processor
784
785 if (returnReduceOr(nPatchFaces))
786 {
787 patches[nValidPatches] = polyPatch::New
788 (
789 patchNames[p],
790 patchDicts[p],
791 nValidPatches,
792 pMesh.boundaryMesh()
793 ).ptr();
794
795 nValidPatches++;
796 }
797 }
798 }
799
800 patches.setSize(nValidPatches);
801
802 pMesh.addPatches(patches);
803
804 return meshPtr;
805}
806
807
808void Foam::conformalVoronoiMesh::checkCellSizing()
809{
810 Info<< "Checking cell sizes..."<< endl;
811
812 timeCheck("Start of Cell Sizing");
813
814 labelList boundaryPts(number_of_finite_cells(), internal);
815 pointField ptsField;
816
817 indexDualVertices(ptsField, boundaryPts);
818
819 // Merge close dual vertices.
820 mergeIdenticalDualVertices(ptsField, boundaryPts);
821
822 autoPtr<polyMesh> meshPtr = createPolyMeshFromPoints(ptsField);
823 const polyMesh& pMesh = meshPtr();
824
825 //pMesh.write();
826
827 // Find cells with poor quality
828 DynamicList<label> checkFaces(identity(pMesh.nFaces()));
829 labelHashSet wrongFaces(pMesh.nFaces()/100);
830
831 Info<< "Running checkMesh on mesh with " << pMesh.nCells()
832 << " cells "<< endl;
833
834 const dictionary& dict
835 = foamyHexMeshControls().foamyHexMeshDict();
836
837 const dictionary& meshQualityDict
838 = dict.subDict("meshQualityControls");
839
840 const scalar maxNonOrtho =
841 meshQualityDict.get<scalar>("maxNonOrtho", keyType::REGEX_RECURSIVE);
842
843 label nWrongFaces = 0;
844
845 if (maxNonOrtho < 180.0 - SMALL)
846 {
848 (
849 false,
850 maxNonOrtho,
851 pMesh,
852 pMesh.cellCentres(),
853 pMesh.faceAreas(),
854 checkFaces,
856 &wrongFaces
857 );
858
859 label nNonOrthogonal = returnReduce(wrongFaces.size(), sumOp<label>());
860
861 Info<< " non-orthogonality > " << maxNonOrtho
862 << " degrees : " << nNonOrthogonal << endl;
863
864 nWrongFaces += nNonOrthogonal;
865 }
866
867 labelHashSet protrudingCells = findOffsetPatchFaces(pMesh, 0.25);
868
869 label nProtrudingCells = protrudingCells.size();
870
871 Info<< " protruding/intruding cells : " << nProtrudingCells << endl;
872
873 nWrongFaces += nProtrudingCells;
874
875// motionSmoother::checkMesh
876// (
877// false,
878// pMesh,
879// meshQualityDict,
880// checkFaces,
881// wrongFaces
882// );
883
884 Info<< " Found total of " << nWrongFaces << " bad faces" << endl;
885
886 {
887 labelHashSet cellsToResizeMap(pMesh.nFaces()/100);
888
889 // Find cells that are attached to the faces in wrongFaces.
890 for (const label facei : wrongFaces)
891 {
892 const label faceOwner = pMesh.faceOwner()[facei];
893 const label faceNeighbour = pMesh.faceNeighbour()[facei];
894
895 cellsToResizeMap.insert(faceOwner);
896 cellsToResizeMap.insert(faceNeighbour);
897 }
898
899 cellsToResizeMap += protrudingCells;
900
901 pointField cellsToResize(cellsToResizeMap.size());
902
903 label count = 0;
904 for (label celli = 0; celli < pMesh.nCells(); ++celli)
905 {
906 if (cellsToResizeMap.found(celli))
907 {
908 cellsToResize[count++] = pMesh.cellCentres()[celli];
909 }
910 }
911
912 Info<< " DISABLED: Automatically re-sizing " << cellsToResize.size()
913 << " cells that are attached to the bad faces: " << endl;
914
915 //cellSizeControl_.setCellSizes(cellsToResize);
916 }
917
918 timeCheck("End of Cell Sizing");
919
920 Info<< "Finished checking cell sizes"<< endl;
921}
922
923
924Foam::labelHashSet Foam::conformalVoronoiMesh::findOffsetPatchFaces
925(
926 const polyMesh& mesh,
927 const scalar allowedOffset
928) const
929{
930 timeCheck("Start findRemainingProtrusionSet");
931
933
934 cellSet offsetBoundaryCells
935 (
936 mesh,
937 "foamyHexMesh_protrudingCells",
938 mesh.nCells()/1000
939 );
940
941 forAll(patches, patchi)
942 {
943 const polyPatch& patch = patches[patchi];
944
945 const faceList& localFaces = patch.localFaces();
946 const pointField& localPoints = patch.localPoints();
947
948 const labelUList& fCell = patch.faceCells();
949
950 forAll(localFaces, pLFI)
951 {
952 const face& f = localFaces[pLFI];
953
954 const Foam::point& faceCentre = f.centre(localPoints);
955
956 const scalar targetSize = targetCellSize(faceCentre);
957
958 pointIndexHit pHit;
959 label surfHit = -1;
960
961 geometryToConformTo_.findSurfaceNearest
962 (
963 faceCentre,
964 sqr(targetSize),
965 pHit,
966 surfHit
967 );
968
969 if
970 (
971 pHit.hit()
972 && (pHit.point().dist(faceCentre) > allowedOffset*targetSize)
973 )
974 {
975 offsetBoundaryCells.insert(fCell[pLFI]);
976 }
977 }
978 }
979
980 if (foamyHexMeshControls().objOutput())
981 {
982 offsetBoundaryCells.write();
983 }
984
985 return offsetBoundaryCells;
986}
987
988
989Foam::labelHashSet Foam::conformalVoronoiMesh::checkPolyMeshQuality
990(
991 const pointField& pts
992) const
993{
994 autoPtr<polyMesh> meshPtr = createPolyMeshFromPoints(pts);
995 polyMesh& pMesh = meshPtr();
996
997 timeCheck("polyMesh created, checking quality");
998
999 labelHashSet wrongFaces(pMesh.nFaces()/100);
1000
1001 DynamicList<label> checkFaces(pMesh.nFaces());
1002
1003 const vectorField& fAreas = pMesh.faceAreas();
1004
1005 scalar faceAreaLimit = SMALL;
1006
1007 forAll(fAreas, fI)
1008 {
1009 if (mag(fAreas[fI]) > faceAreaLimit)
1010 {
1011 checkFaces.append(fI);
1012 }
1013 }
1014
1015 Info<< nl << "Excluding "
1016 << returnReduce(fAreas.size() - checkFaces.size(), sumOp<label>())
1017 << " faces from check, < " << faceAreaLimit << " area" << endl;
1018
1019 const dictionary& dict
1020 = foamyHexMeshControls().foamyHexMeshDict();
1021
1022 const dictionary& meshQualityDict
1023 = dict.subDict("meshQualityControls");
1024
1026 (
1027 false,
1028 pMesh,
1029 meshQualityDict,
1030 checkFaces,
1031 wrongFaces
1032 );
1033
1034 {
1035 // Check for cells with more than 1 but fewer than 4 faces
1036 label nInvalidPolyhedra = 0;
1037
1038 const cellList& cells = pMesh.cells();
1039
1040 forAll(cells, cI)
1041 {
1042 if (cells[cI].size() < 4 && cells[cI].size() > 0)
1043 {
1044 // Pout<< "cell " << cI << " " << cells[cI]
1045 // << " has " << cells[cI].size() << " faces."
1046 // << endl;
1047
1048 nInvalidPolyhedra++;
1049
1050 wrongFaces.insert(cells[cI]);
1051 }
1052 }
1053
1054 Info<< " cells with more than 1 but fewer than 4 faces : "
1055 << returnReduce(nInvalidPolyhedra, sumOp<label>())
1056 << endl;
1057
1058 // Check for cells with one internal face only
1059
1060 labelList nInternalFaces(pMesh.nCells(), Zero);
1061
1062 for (label fI = 0; fI < pMesh.nInternalFaces(); fI++)
1063 {
1064 nInternalFaces[pMesh.faceOwner()[fI]]++;
1065 nInternalFaces[pMesh.faceNeighbour()[fI]]++;
1066 }
1067
1068 const polyBoundaryMesh& patches = pMesh.boundaryMesh();
1069
1070 forAll(patches, patchi)
1071 {
1072 if (patches[patchi].coupled())
1073 {
1074 const labelUList& owners = patches[patchi].faceCells();
1075
1076 forAll(owners, i)
1077 {
1078 nInternalFaces[owners[i]]++;
1079 }
1080 }
1081 }
1082
1083 label oneInternalFaceCells = 0;
1084
1085 forAll(nInternalFaces, cI)
1086 {
1087 if (nInternalFaces[cI] <= 1)
1088 {
1089 oneInternalFaceCells++;
1090 wrongFaces.insert(cells[cI]);
1091 }
1092 }
1093
1094 Info<< " cells with with zero or one non-boundary face : "
1095 << returnReduce(oneInternalFaceCells, sumOp<label>())
1096 << endl;
1097 }
1098
1099
1100 bitSet ptToBeLimited(pts.size(), false);
1101
1102 for (const label facei : wrongFaces)
1103 {
1104 const face f = pMesh.faces()[facei];
1105
1106 ptToBeLimited.set(f);
1107 }
1108
1109 // // Limit connected cells
1110
1111 // labelHashSet limitCells(pMesh.nCells()/100);
1112
1113 // const labelListList& ptCells = pMesh.pointCells();
1114
1115 // for (const label facei : wrongFaces)
1116 // {
1117 // const face f = pMesh.faces()[facei];
1118
1119 // forAll(f, fPtI)
1120 // {
1121 // label ptI = f[fPtI];
1122 // const labelList& pC = ptCells[ptI];
1123 // limitCells.insert(pC);
1124 // }
1125 // }
1126
1127 // const labelListList& cellPts = pMesh.cellPoints();
1128
1129 // for (const label celli : limitCells)
1130 // {
1131 // const labelList& cP = cellPts[celli];
1132
1133 // ptToBeLimited.set(cP);
1134 // }
1135
1136
1137 // Apply Delaunay cell filterCounts and determine the maximum
1138 // overall filterCount
1139
1140 label maxFilterCount = 0;
1141
1142 for
1143 (
1144 Delaunay::Finite_cells_iterator cit = finite_cells_begin();
1145 cit != finite_cells_end();
1146 ++cit
1147 )
1148 {
1149 label cI = cit->cellIndex();
1150
1151 if (cI >= 0)
1152 {
1153 if (ptToBeLimited.test(cI))
1154 {
1155 cit->filterCount()++;
1156 }
1157
1158 if (cit->filterCount() > maxFilterCount)
1159 {
1160 maxFilterCount = cit->filterCount();
1161 }
1162 }
1163 }
1164
1165 Info<< nl << "Maximum number of filter limits applied: "
1166 << returnReduce(maxFilterCount, maxOp<label>()) << endl;
1167
1168 return wrongFaces;
1169}
1170
1171
1172Foam::label Foam::conformalVoronoiMesh::classifyBoundaryPoint
1173(
1174 Cell_handle cit
1175) const
1176{
1177 if (cit->boundaryDualVertex())
1178 {
1179 if (cit->featurePointDualVertex())
1180 {
1181 return featurePoint;
1182 }
1183 else if (cit->featureEdgeDualVertex())
1184 {
1185 return featureEdge;
1186 }
1187 else
1188 {
1189 return surface;
1190 }
1191 }
1192 else if (cit->baffleSurfaceDualVertex())
1193 {
1194 return surface;
1195 }
1196 else if (cit->baffleEdgeDualVertex())
1197 {
1198 return featureEdge;
1199 }
1200 else
1201 {
1202 return internal;
1203 }
1204}
1205
1206
1207void Foam::conformalVoronoiMesh::indexDualVertices
1208(
1209 pointField& pts,
1210 labelList& boundaryPts
1211)
1212{
1213 // Indexing Delaunay cells, which are the dual vertices
1214
1215 this->resetCellCount();
1216
1217 label nConstrainedVertices = 0;
1218 if (foamyHexMeshControls().guardFeaturePoints())
1219 {
1220 for
1221 (
1222 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1223 vit != finite_vertices_end();
1224 ++vit
1225 )
1226 {
1227 if (vit->constrained())
1228 {
1229 vit->index() = number_of_finite_cells() + nConstrainedVertices;
1230 nConstrainedVertices++;
1231 }
1232 }
1233 }
1234
1235 pts.setSize(number_of_finite_cells() + nConstrainedVertices);
1236 boundaryPts.setSize
1237 (
1238 number_of_finite_cells() + nConstrainedVertices,
1239 internal
1240 );
1241
1242 if (foamyHexMeshControls().guardFeaturePoints())
1243 {
1244 nConstrainedVertices = 0;
1245 for
1246 (
1247 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1248 vit != finite_vertices_end();
1249 ++vit
1250 )
1251 {
1252 if (vit->constrained())
1253 {
1254 pts[number_of_finite_cells() + nConstrainedVertices] =
1255 topoint(vit->point());
1256
1257 boundaryPts[number_of_finite_cells() + nConstrainedVertices] =
1258 constrained;
1259
1260 nConstrainedVertices++;
1261 }
1262 }
1263 }
1264
1265 //OBJstream snapping1("snapToSurface1.obj");
1266 //OBJstream snapping2("snapToSurface2.obj");
1267 //OFstream tetToSnapTo("tetsToSnapTo.obj");
1268
1269 for
1270 (
1271 Delaunay::Finite_cells_iterator cit = finite_cells_begin();
1272 cit != finite_cells_end();
1273 ++cit
1274 )
1275 {
1276// if (tetrahedron(cit).volume() == 0)
1277// {
1278// Pout<< "ZERO VOLUME TET" << endl;
1279// Pout<< cit->info();
1280// Pout<< "Dual = " << cit->dual();
1281// }
1282
1283 if (!cit->hasFarPoint())
1284 {
1285 cit->cellIndex() = getNewCellIndex();
1286
1287 // For nearly coplanar Delaunay cells that are present on different
1288 // processors the result of the circumcentre calculation depends on
1289 // the ordering of the vertices, so synchronise it across processors
1290
1291 if (Pstream::parRun() && cit->parallelDualVertex())
1292 {
1293 typedef CGAL::Exact_predicates_exact_constructions_kernel Exact;
1294 typedef CGAL::Point_3<Exact> ExactPoint;
1295
1296 List<labelPair> cellVerticesPair(4);
1297 List<ExactPoint> cellVertices(4);
1298
1299 for (label vI = 0; vI < 4; ++vI)
1300 {
1301 cellVerticesPair[vI] = labelPair
1302 (
1303 cit->vertex(vI)->procIndex(),
1304 cit->vertex(vI)->index()
1305 );
1306
1307 cellVertices[vI] = ExactPoint
1308 (
1309 cit->vertex(vI)->point().x(),
1310 cit->vertex(vI)->point().y(),
1311 cit->vertex(vI)->point().z()
1312 );
1313 }
1314
1315 // Sort the vertices so that they will be in the same order on
1316 // each processor
1317 labelList oldToNew(sortedOrder(cellVerticesPair));
1318 oldToNew = invert(oldToNew.size(), oldToNew);
1319 inplaceReorder(oldToNew, cellVertices);
1320
1321 ExactPoint synchronisedDual = CGAL::circumcenter
1322 (
1323 cellVertices[0],
1324 cellVertices[1],
1325 cellVertices[2],
1326 cellVertices[3]
1327 );
1328
1329 pts[cit->cellIndex()] = Foam::point
1330 (
1331 CGAL::to_double(synchronisedDual.x()),
1332 CGAL::to_double(synchronisedDual.y()),
1333 CGAL::to_double(synchronisedDual.z())
1334 );
1335 }
1336 else
1337 {
1338 pts[cit->cellIndex()] = cit->dual();
1339 }
1340
1341 // Feature point snapping
1342 if (foamyHexMeshControls().snapFeaturePoints())
1343 {
1344 if (cit->featurePointDualVertex())
1345 {
1346 pointFromPoint dual = cit->dual();
1347
1348 pointIndexHit fpHit;
1349 label featureHit;
1350
1351 // Find nearest feature point and compare
1352 geometryToConformTo_.findFeaturePointNearest
1353 (
1354 dual,
1355 sqr(targetCellSize(dual)),
1356 fpHit,
1357 featureHit
1358 );
1359
1360 if (fpHit.hit())
1361 {
1362 if (debug)
1363 {
1364 Info<< "Dual = " << dual << nl
1365 << " Nearest = " << fpHit.point() << endl;
1366 }
1367
1368 pts[cit->cellIndex()] = fpHit.point();
1369 }
1370 }
1371 }
1372
1373// {
1374// // Snapping points far outside
1375// if (cit->boundaryDualVertex() && !cit->parallelDualVertex())
1376// {
1377// pointFromPoint dual = cit->dual();
1378//
1379// pointIndexHit hitInfo;
1380// label surfHit;
1381//
1382// // Find nearest surface point
1383// geometryToConformTo_.findSurfaceNearest
1384// (
1385// dual,
1386// sqr(targetCellSize(dual)),
1387// hitInfo,
1388// surfHit
1389// );
1390//
1391// if (!hitInfo.hit())
1392// {
1393// // Project dual to nearest point on tet
1394//
1395// tetPointRef tet
1396// (
1397// topoint(cit->vertex(0)->point()),
1398// topoint(cit->vertex(1)->point()),
1399// topoint(cit->vertex(2)->point()),
1400// topoint(cit->vertex(3)->point())
1401// );
1402//
1403// pointFromPoint nearestPointOnTet =
1404// tet.nearestPoint(dual).point();
1405//
1406// // Get nearest point on surface from tet.
1407// geometryToConformTo_.findSurfaceNearest
1408// (
1409// nearestPointOnTet,
1410// sqr(targetCellSize(nearestPointOnTet)),
1411// hitInfo,
1412// surfHit
1413// );
1414//
1415// vector snapDir = nearestPointOnTet - dual;
1416// snapDir /= mag(snapDir) + SMALL;
1417//
1418// drawDelaunayCell(tetToSnapTo, cit, offset);
1419// offset += 1;
1420//
1421// vectorField norm(1);
1422// allGeometry_[surfHit].getNormal
1423// (
1424// List<pointIndexHit>(1, hitInfo),
1425// norm
1426// );
1427// norm[0] /= mag(norm[0]) + SMALL;
1428//
1429// if
1430// (
1431// hitInfo.hit()
1432// && (mag(snapDir & norm[0]) > 0.5)
1433// )
1434// {
1435// snapping1.writeLine
1436// (
1437// dual,
1438// nearestPointOnTet
1439// );
1440// snapping2.writeLine
1441// (
1442// nearestPointOnTet,
1443// hitInfo.point()
1444// );
1445//
1446// pts[cit->cellIndex()] = hitInfo.point();
1447// }
1448// }
1449// }
1450// }
1451
1452 boundaryPts[cit->cellIndex()] = classifyBoundaryPoint(cit);
1453 }
1454 else
1455 {
1456 cit->cellIndex() = Cb::ctFar;
1457 }
1458 }
1459
1460 //pts.setSize(this->cellCount());
1461
1462 //boundaryPts.setSize(this->cellCount());
1463}
1464
1465
1466void Foam::conformalVoronoiMesh::reindexDualVertices
1467(
1468 const Map<label>& dualPtIndexMap,
1469 labelList& boundaryPts
1470)
1471{
1472 for
1473 (
1474 Delaunay::Finite_cells_iterator cit = finite_cells_begin();
1475 cit != finite_cells_end();
1476 ++cit
1477 )
1478 {
1479 if (dualPtIndexMap.found(cit->cellIndex()))
1480 {
1481 cit->cellIndex() = dualPtIndexMap[cit->cellIndex()];
1482 boundaryPts[cit->cellIndex()] =
1483 max
1484 (
1485 boundaryPts[cit->cellIndex()],
1486 boundaryPts[dualPtIndexMap[cit->cellIndex()]]
1487 );
1488 }
1489 }
1490}
1491
1492
1493Foam::label Foam::conformalVoronoiMesh::createPatchInfo
1494(
1497) const
1498{
1499 patchNames = geometryToConformTo_.patchNames();
1500
1502
1503 const PtrList<dictionary>& patchInfo = geometryToConformTo_.patchInfo();
1504
1505 forAll(patchNames, patchi)
1506 {
1507 if (patchInfo.set(patchi))
1508 {
1509 patchDicts.set(patchi, new dictionary(patchInfo[patchi]));
1510 }
1511 else
1512 {
1513 patchDicts.set(patchi, new dictionary());
1514 patchDicts[patchi].set
1515 (
1516 "type",
1517 wallPolyPatch::typeName
1518 );
1519 }
1520 }
1521
1523 label defaultPatchIndex = patchNames.size() - 1;
1524 patchNames[defaultPatchIndex] = "foamyHexMesh_defaultPatch";
1525 patchDicts.set(defaultPatchIndex, new dictionary());
1526 patchDicts[defaultPatchIndex].set
1527 (
1528 "type",
1529 wallPolyPatch::typeName
1530 );
1531
1532 label nProcPatches = 0;
1533
1534 if (Pstream::parRun())
1535 {
1536 List<boolList> procUsedList
1537 (
1539 boolList(Pstream::nProcs(), false)
1540 );
1541
1542 boolList& procUsed = procUsedList[Pstream::myProcNo()];
1543
1544 // Determine which processor patches are required
1545 for
1546 (
1547 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1548 vit != finite_vertices_end();
1549 vit++
1550 )
1551 {
1552 // This test is not sufficient if one of the processors does
1553 // not receive a referred vertex from another processor, but does
1554 // send one to the other processor.
1555 if (vit->referred())
1556 {
1557 procUsed[vit->procIndex()] = true;
1558 }
1559 }
1560
1561 // Because the previous test was insufficient, combine the lists.
1562 Pstream::allGatherList(procUsedList);
1563
1564 forAll(procUsedList, proci)
1565 {
1566 if (proci != Pstream::myProcNo())
1567 {
1568 if (procUsedList[proci][Pstream::myProcNo()])
1569 {
1570 procUsed[proci] = true;
1571 }
1572 }
1573 }
1574
1575 forAll(procUsed, pUI)
1576 {
1577 if (procUsed[pUI])
1578 {
1579 nProcPatches++;
1580 }
1581 }
1582
1583 label nNonProcPatches = patchNames.size();
1584 label nTotalPatches = nNonProcPatches + nProcPatches;
1585
1586 patchNames.setSize(nTotalPatches);
1587 patchDicts.setSize(nTotalPatches);
1588 for (label pI = nNonProcPatches; pI < nTotalPatches; ++pI)
1589 {
1590 patchDicts.set(pI, new dictionary());
1591 }
1592
1593 label procAddI = 0;
1594
1595 forAll(procUsed, pUI)
1596 {
1597 if (procUsed[pUI])
1598 {
1599 patchNames[nNonProcPatches + procAddI] =
1601
1602 patchDicts[nNonProcPatches + procAddI].set
1603 (
1604 "type",
1605 processorPolyPatch::typeName
1606 );
1607
1608 patchDicts[nNonProcPatches + procAddI].set
1609 (
1610 "myProcNo",
1612 );
1613
1614 patchDicts[nNonProcPatches + procAddI].set("neighbProcNo", pUI);
1615
1616 procAddI++;
1617 }
1618 }
1619 }
1620
1621 return defaultPatchIndex;
1622}
1623
1624
1625Foam::vector Foam::conformalVoronoiMesh::calcSharedPatchNormal
1626(
1627 Cell_handle c1,
1628 Cell_handle c2
1629) const
1630{
1631 List<Foam::point> patchEdge(2, point::max);
1632
1633 // Get shared Facet
1634 for (label cI = 0; cI < 4; ++cI)
1635 {
1636 if (c1->neighbor(cI) != c2 && !c1->vertex(cI)->constrained())
1637 {
1638 if (c1->vertex(cI)->internalBoundaryPoint())
1639 {
1640 patchEdge[0] = topoint(c1->vertex(cI)->point());
1641 }
1642 else
1643 {
1644 patchEdge[1] = topoint(c1->vertex(cI)->point());
1645 }
1646 }
1647 }
1648
1649 Info<< " " << patchEdge << endl;
1650
1651 return vector(patchEdge[1] - patchEdge[0]);
1652}
1653
1654
1655bool Foam::conformalVoronoiMesh::boundaryDualFace
1656(
1657 Cell_handle c1,
1658 Cell_handle c2
1659) const
1660{
1661 label nInternal = 0;
1662 label nExternal = 0;
1663
1664 for (label cI = 0; cI < 4; ++cI)
1665 {
1666 if (c1->neighbor(cI) != c2 && !c1->vertex(cI)->constrained())
1667 {
1668 if (c1->vertex(cI)->internalBoundaryPoint())
1669 {
1670 nInternal++;
1671 }
1672 else if (c1->vertex(cI)->externalBoundaryPoint())
1673 {
1674 nExternal++;
1675 }
1676 }
1677 }
1678
1679 Info<< "in = " << nInternal << " out = " << nExternal << endl;
1680
1681 return (nInternal == 1 && nExternal == 1);
1682}
1683
1684
1685void Foam::conformalVoronoiMesh::createFacesOwnerNeighbourAndPatches
1686(
1687 const pointField& pts,
1688 faceList& faces,
1689 labelList& owner,
1690 labelList& neighbour,
1693 labelListList& patchPointPairSlaves,
1694 bitSet& boundaryFacesToRemove,
1695 bool includeEmptyPatches
1696) const
1697{
1698 const label defaultPatchIndex = createPatchInfo(patchNames, patchDicts);
1699
1700 const label nPatches = patchNames.size();
1701
1702 labelList procNeighbours(nPatches);
1703 forAll(procNeighbours, patchi)
1704 {
1705 procNeighbours[patchi] =
1706 patchDicts[patchi].getOrDefault<label>("neighbProcNo", -1);
1707 }
1708
1709 List<DynamicList<face>> patchFaces(nPatches);
1710 List<DynamicList<label>> patchOwners(nPatches);
1711 // Per patch face the index of the slave node of the point pair
1712 List<DynamicList<label>> patchPPSlaves(nPatches);
1713 List<DynamicList<bool>> indirectPatchFace(nPatches);
1714
1715
1716 faces.setSize(number_of_finite_edges());
1717 owner.setSize(number_of_finite_edges());
1718 neighbour.setSize(number_of_finite_edges());
1719 boundaryFacesToRemove.setSize(number_of_finite_edges(), false);
1720
1721 labelPairPairDynListList procPatchSortingIndex(nPatches);
1722
1723 label dualFacei = 0;
1724
1725 if (foamyHexMeshControls().guardFeaturePoints())
1726 {
1727 OBJstream startCellStr("startingCell.obj");
1728 OBJstream featurePointFacesStr("ftPtFaces.obj");
1729 OBJstream featurePointDualsStr("ftPtDuals.obj");
1730 OFstream cellStr("vertexCells.obj");
1731
1732 label vcount = 1;
1733
1734 for
1735 (
1736 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1737 vit != finite_vertices_end();
1738 ++vit
1739 )
1740 {
1741 if (vit->constrained())
1742 {
1743 // Find a starting cell
1744 std::list<Cell_handle> vertexCells;
1745 finite_incident_cells(vit, std::back_inserter(vertexCells));
1746
1747 Cell_handle startCell;
1748
1749 for
1750 (
1751 std::list<Cell_handle>::iterator vcit = vertexCells.begin();
1752 vcit != vertexCells.end();
1753 ++vcit
1754 )
1755 {
1756 if ((*vcit)->featurePointExternalCell())
1757 {
1758 startCell = *vcit;
1759 }
1760
1761 if ((*vcit)->real())
1762 {
1763 featurePointDualsStr.writeLine
1764 (
1765 topoint(vit->point()),
1766 (*vcit)->dual()
1767 );
1768 }
1769 }
1770
1771 // Error if startCell is null
1772 if (startCell == nullptr)
1773 {
1774 Pout<< "Start cell is null!" << endl;
1775 }
1776
1777 // Need to pick a direction to walk in
1778 Cell_handle vc1 = startCell;
1779 Cell_handle vc2;
1780
1781 Info<< "c1 index = " << vc1->cellIndex() << " "
1782 << vc1->dual() << endl;
1783
1784 for (label cI = 0; cI < 4; ++cI)
1785 {
1786 Info<< "c1 = " << cI << " "
1787 << vc1->neighbor(cI)->cellIndex() << " v = "
1788 << vc1->neighbor(cI)->dual() << endl;
1789
1790 Info<< vc1->vertex(cI)->info();
1791 }
1792
1793 Cell_handle nextCell;
1794
1795 for (label cI = 0; cI < 4; ++cI)
1796 {
1797 if (vc1->vertex(cI)->externalBoundaryPoint())
1798 {
1799 vc2 = vc1->neighbor(cI);
1800
1801 Info<< " c2 is neighbor "
1802 << vc2->cellIndex()
1803 << " of c1" << endl;
1804
1805 for (label cI = 0; cI < 4; ++cI)
1806 {
1807 Info<< " c2 = " << cI << " "
1808 << vc2->neighbor(cI)->cellIndex() << " v = "
1809 << vc2->vertex(cI)->index() << endl;
1810 }
1811
1812 face f(3);
1813 f[0] = vit->index();
1814 f[1] = vc1->cellIndex();
1815 f[2] = vc2->cellIndex();
1816
1817 Info<< "f " << f << endl;
1818 forAll(f, pI)
1819 {
1820 Info<< " " << pts[f[pI]] << endl;
1821 }
1822
1823 vector correctNormal = calcSharedPatchNormal(vc1, vc2);
1824 correctNormal.normalise();
1825
1826 Info<< " cN " << correctNormal << endl;
1827
1828 vector fN = f.areaNormal(pts);
1829
1830 if (mag(fN) < SMALL)
1831 {
1832 nextCell = vc2;
1833 continue;
1834 }
1835
1836 fN.normalise();
1837 Info<< " fN " << fN << endl;
1838
1839 if ((fN & correctNormal) > 0)
1840 {
1841 nextCell = vc2;
1842 break;
1843 }
1844 }
1845 }
1846
1847 vc2 = nextCell;
1848
1849 label own = vit->index();
1850 face f(3);
1851 f[0] = own;
1852
1853 Info<< "Start walk from " << vc1->cellIndex()
1854 << " to " << vc2->cellIndex() << endl;
1855
1856 // Walk while not at start cell
1857
1858 label iter = 0;
1859 do
1860 {
1861 Info<< " Walk from " << vc1->cellIndex()
1862 << " " << vc1->dual()
1863 << " to " << vc2->cellIndex()
1864 << " " << vc2->dual()
1865 << endl;
1866
1867 startCellStr.writeLine(vc1->dual(), vc2->dual());
1868
1869 // Get patch by getting face between cells and the two
1870 // points on the face that are not the feature vertex
1871 label patchIndex =
1872 geometryToConformTo_.findPatch
1873 (
1874 topoint(vit->point())
1875 );
1876
1877 f[1] = vc1->cellIndex();
1878 f[2] = vc2->cellIndex();
1879
1880 patchFaces[patchIndex].append(f);
1881 patchOwners[patchIndex].append(own);
1882 patchPPSlaves[patchIndex].append(own);
1883
1884 // Find next cell
1885 Cell_handle nextCell;
1886
1887 Info<< " c1 vertices " << vc2->dual() << endl;
1888 for (label cI = 0; cI < 4; ++cI)
1889 {
1890 Info<< " " << vc2->vertex(cI)->info();
1891 }
1892 Info<< " c1 neighbour vertices " << endl;
1893 for (label cI = 0; cI < 4; ++cI)
1894 {
1895 if
1896 (
1897 !vc2->vertex(cI)->constrained()
1898 && vc2->neighbor(cI) != vc1
1899 && !is_infinite(vc2->neighbor(cI))
1900 &&
1901 (
1902 vc2->neighbor(cI)->featurePointExternalCell()
1903 || vc2->neighbor(cI)->featurePointInternalCell()
1904 )
1905 && vc2->neighbor(cI)->hasConstrainedPoint()
1906 )
1907 {
1909 (
1910 cellStr,
1911 vc2->neighbor(cI),
1912 vcount++
1913 );
1914
1915 Info<< " neighbour " << cI << " "
1916 << vc2->neighbor(cI)->dual() << endl;
1917 for (label I = 0; I < 4; ++I)
1918 {
1919 Info<< " "
1920 << vc2->neighbor(cI)->vertex(I)->info();
1921 }
1922 }
1923 }
1924
1925 for (label cI = 0; cI < 4; ++cI)
1926 {
1927 if
1928 (
1929 !vc2->vertex(cI)->constrained()
1930 && vc2->neighbor(cI) != vc1
1931 && !is_infinite(vc2->neighbor(cI))
1932 &&
1933 (
1934 vc2->neighbor(cI)->featurePointExternalCell()
1935 || vc2->neighbor(cI)->featurePointInternalCell()
1936 )
1937 && vc2->neighbor(cI)->hasConstrainedPoint()
1938 )
1939 {
1940 // check if shared edge is internal/internal
1941 if (boundaryDualFace(vc2, vc2->neighbor(cI)))
1942 {
1943 nextCell = vc2->neighbor(cI);
1944 break;
1945 }
1946 }
1947 }
1948
1949 vc1 = vc2;
1950 vc2 = nextCell;
1951
1952 iter++;
1953 } while (vc1 != startCell && iter < 100);
1954 }
1955 }
1956 }
1957
1958 for
1959 (
1960 Delaunay::Finite_edges_iterator eit = finite_edges_begin();
1961 eit != finite_edges_end();
1962 ++eit
1963 )
1964 {
1965 Cell_handle c = eit->first;
1966 Vertex_handle vA = c->vertex(eit->second);
1967 Vertex_handle vB = c->vertex(eit->third);
1968
1969 if (vA->constrained() && vB->constrained())
1970 {
1971 continue;
1972 }
1973
1974 if
1975 (
1976 (vA->constrained() && vB->internalOrBoundaryPoint())
1977 || (vB->constrained() && vA->internalOrBoundaryPoint())
1978 )
1979 {
1980 face newDualFace = buildDualFace(eit);
1981
1982 label own = -1;
1983 label nei = -1;
1984
1985 if (ownerAndNeighbour(vA, vB, own, nei))
1986 {
1987 reverse(newDualFace);
1988 }
1989
1990 // internal face
1991 faces[dualFacei] = newDualFace;
1992 owner[dualFacei] = own;
1993 neighbour[dualFacei] = nei;
1994
1995 dualFacei++;
1996 }
1997 else if
1998 (
1999 (vA->internalOrBoundaryPoint() && !vA->referred())
2000 || (vB->internalOrBoundaryPoint() && !vB->referred())
2001 )
2002 {
2003 if
2004 (
2005 (vA->internalPoint() && vB->externalBoundaryPoint())
2006 || (vB->internalPoint() && vA->externalBoundaryPoint())
2007 )
2008 {
2009 Cell_circulator ccStart = incident_cells(*eit);
2010 Cell_circulator cc1 = ccStart;
2011 Cell_circulator cc2 = cc1;
2012
2013 cc2++;
2014
2015 bool skipEdge = false;
2016
2017 do
2018 {
2019 if
2020 (
2021 cc1->hasFarPoint() || cc2->hasFarPoint()
2022 || is_infinite(cc1) || is_infinite(cc2)
2023 )
2024 {
2025 Pout<< "Ignoring edge between internal and external: "
2026 << vA->info()
2027 << vB->info();
2028
2029 skipEdge = true;
2030 break;
2031 }
2032
2033 cc1++;
2034 cc2++;
2035
2036 } while (cc1 != ccStart);
2037
2038
2039 // Do not create faces if the internal point is outside!
2040 // This occurs because the internal point is not determined to
2041 // be outside in the inside/outside test. This is most likely
2042 // due to the triangle.nearestPointClassify test not returning
2043 // edge/point as the nearest type.
2044
2045 if (skipEdge)
2046 {
2047 continue;
2048 }
2049 }
2050
2051 face newDualFace = buildDualFace(eit);
2052
2053 if (newDualFace.size() >= 3)
2054 {
2055 label own = -1;
2056 label nei = -1;
2057
2058 if (ownerAndNeighbour(vA, vB, own, nei))
2059 {
2060 reverse(newDualFace);
2061 }
2062
2063 label patchIndex = -1;
2064
2065 pointFromPoint ptA = topoint(vA->point());
2066 pointFromPoint ptB = topoint(vB->point());
2067
2068 if (nei == -1)
2069 {
2070 // boundary face
2071
2072 if (isProcBoundaryEdge(eit))
2073 {
2074 // One (and only one) of the points is an internal
2075 // point from another processor
2076
2077 label procIndex = max(vA->procIndex(), vB->procIndex());
2078
2079 patchIndex = max
2080 (
2081 procNeighbours.find(vA->procIndex()),
2082 procNeighbours.find(vB->procIndex())
2083 );
2084
2085 // The lower processor index is the owner of the
2086 // two for the purpose of sorting the patch faces.
2087
2088 if (Pstream::myProcNo() < procIndex)
2089 {
2090 // Use this processor's vertex index as the master
2091 // for sorting
2092
2093 DynamicList<labelPairPair>& sortingIndex =
2094 procPatchSortingIndex[patchIndex];
2095
2096 if (vB->internalOrBoundaryPoint() && vB->referred())
2097 {
2098 sortingIndex.append
2099 (
2101 (
2102 labelPair(vA->index(), vA->procIndex()),
2103 labelPair(vB->index(), vB->procIndex())
2104 )
2105 );
2106 }
2107 else
2108 {
2109 sortingIndex.append
2110 (
2112 (
2113 labelPair(vB->index(), vB->procIndex()),
2114 labelPair(vA->index(), vA->procIndex())
2115 )
2116 );
2117 }
2118 }
2119 else
2120 {
2121 // Use the other processor's vertex index as the
2122 // master for sorting
2123
2124 DynamicList<labelPairPair>& sortingIndex =
2125 procPatchSortingIndex[patchIndex];
2126
2127 if (vA->internalOrBoundaryPoint() && vA->referred())
2128 {
2129 sortingIndex.append
2130 (
2132 (
2133 labelPair(vA->index(), vA->procIndex()),
2134 labelPair(vB->index(), vB->procIndex())
2135 )
2136 );
2137 }
2138 else
2139 {
2140 sortingIndex.append
2141 (
2143 (
2144 labelPair(vB->index(), vB->procIndex()),
2145 labelPair(vA->index(), vA->procIndex())
2146 )
2147 );
2148 }
2149 }
2150
2151// Pout<< ptA << " " << ptB
2152// << " proc indices "
2153// << vA->procIndex() << " " << vB->procIndex()
2154// << " indices " << vA->index()
2155// << " " << vB->index()
2156// << " my proc " << Pstream::myProcNo()
2157// << " addedIndex "
2158// << procPatchSortingIndex[patchIndex].last()
2159// << endl;
2160 }
2161 else
2162 {
2163 patchIndex = geometryToConformTo_.findPatch(ptA, ptB);
2164 }
2165
2166 if (patchIndex == -1)
2167 {
2168 // Did not find a surface patch between
2169 // between Dv pair, finding nearest patch
2170
2171// Pout<< "Did not find a surface patch between "
2172// << "for face, finding nearest patch to"
2173// << 0.5*(ptA + ptB) << endl;
2174
2175 patchIndex = geometryToConformTo_.findPatch
2176 (
2177 0.5*(ptA + ptB)
2178 );
2179 }
2180
2181 patchFaces[patchIndex].append(newDualFace);
2182 patchOwners[patchIndex].append(own);
2183
2184 // If the two vertices are a pair, then the patch face is
2185 // a desired one.
2186 if
2187 (
2188 vA->boundaryPoint() && vB->boundaryPoint()
2189 && !ptPairs_.isPointPair(vA, vB)
2190 && !ftPtConformer_.featurePointPairs().isPointPair(vA, vB)
2191 )
2192 {
2193 indirectPatchFace[patchIndex].append(true);
2194 }
2195 else
2196 {
2197 indirectPatchFace[patchIndex].append(false);
2198 }
2199
2200 // Store the non-internal or boundary point
2201 if (vA->internalOrBoundaryPoint())
2202 {
2203 patchPPSlaves[patchIndex].append(vB->index());
2204 }
2205 else
2206 {
2207 patchPPSlaves[patchIndex].append(vA->index());
2208 }
2209 }
2210 else
2211 {
2212 if
2213 (
2214 !vA->boundaryPoint()
2215 || !vB->boundaryPoint()
2216 || ptPairs_.isPointPair(vA, vB)
2217 || ftPtConformer_.featurePointPairs().isPointPair(vA, vB)
2218 )
2219 {
2220 patchIndex = geometryToConformTo_.findPatch(ptA, ptB);
2221 }
2222
2223 if
2224 (
2225 patchIndex != -1
2226 && geometryToConformTo_.patchInfo().set(patchIndex)
2227 )
2228 {
2229 // baffle faces
2230
2231 patchFaces[patchIndex].append(newDualFace);
2232 patchOwners[patchIndex].append(own);
2233 indirectPatchFace[patchIndex].append(false);
2234
2235 reverse(newDualFace);
2236
2237 patchFaces[patchIndex].append(newDualFace);
2238 patchOwners[patchIndex].append(nei);
2239 indirectPatchFace[patchIndex].append(false);
2240
2241 if
2242 (
2243 labelPair(vB->index(), vB->procIndex())
2244 < labelPair(vA->index(), vA->procIndex())
2245 )
2246 {
2247 patchPPSlaves[patchIndex].append(vB->index());
2248 patchPPSlaves[patchIndex].append(vB->index());
2249 }
2250 else
2251 {
2252 patchPPSlaves[patchIndex].append(vA->index());
2253 patchPPSlaves[patchIndex].append(vA->index());
2254 }
2255
2256 }
2257 else
2258 {
2259 // internal face
2260 faces[dualFacei] = newDualFace;
2261 owner[dualFacei] = own;
2262 neighbour[dualFacei] = nei;
2263
2264 dualFacei++;
2265 }
2266 }
2267 }
2268 }
2269 }
2270
2271 if (!patchFaces[defaultPatchIndex].empty())
2272 {
2273 Pout<< nl << patchFaces[defaultPatchIndex].size()
2274 << " faces were not able to have their patch determined from "
2275 << "the surface. "
2276 << nl << "Adding to patch " << patchNames[defaultPatchIndex]
2277 << endl;
2278 }
2279
2280 label nInternalFaces = dualFacei;
2281
2282 faces.setSize(nInternalFaces);
2283 owner.setSize(nInternalFaces);
2284 neighbour.setSize(nInternalFaces);
2285
2286 timeCheck("polyMesh quality checked");
2287
2288 sortFaces(faces, owner, neighbour);
2289
2290 sortProcPatches
2291 (
2292 patchFaces,
2293 patchOwners,
2294 patchPPSlaves,
2295 procPatchSortingIndex
2296 );
2297
2298 timeCheck("faces, owner, neighbour sorted");
2299
2300 addPatches
2301 (
2302 nInternalFaces,
2303 faces,
2304 owner,
2305 patchDicts,
2306 boundaryFacesToRemove,
2307 patchFaces,
2308 patchOwners,
2309 indirectPatchFace
2310 );
2311
2312 // Return patchPointPairSlaves.setSize(nPatches);
2313 patchPointPairSlaves.setSize(nPatches);
2314 forAll(patchPPSlaves, patchi)
2315 {
2316 patchPointPairSlaves[patchi].transfer(patchPPSlaves[patchi]);
2317 }
2318
2319 if (foamyHexMeshControls().objOutput())
2320 {
2321 Info<< "Writing processor interfaces" << endl;
2322
2323 forAll(patchDicts, nbI)
2324 {
2325 if (patchFaces[nbI].size() > 0)
2326 {
2327 const label neighbour =
2328 patchDicts[nbI].getOrDefault<label>("neighbProcNo", -1);
2329
2330 faceList procPatchFaces = patchFaces[nbI];
2331
2332 // Reverse faces as it makes it easier to analyse the output
2333 // using a diff
2334 if (neighbour < Pstream::myProcNo())
2335 {
2336 forAll(procPatchFaces, fI)
2337 {
2338 procPatchFaces[fI].flip();
2339 }
2340 }
2341
2342 if (neighbour != -1)
2343 {
2344 word fName =
2345 "processor_"
2347 + "_to_"
2348 + name(neighbour)
2349 + "_interface.obj";
2350
2352 (
2353 time().path()/fName,
2354 *this,
2355 procPatchFaces
2356 );
2357 }
2358 }
2359 }
2360 }
2361}
2362
2363
2364void Foam::conformalVoronoiMesh::sortFaces
2365(
2366 faceList& faces,
2367 labelList& owner,
2368 labelList& neighbour
2369) const
2370{
2371 // Upper triangular order:
2372 // + owner is sorted in ascending cell order
2373 // + within each block of equal value for owner, neighbour is sorted in
2374 // ascending cell order.
2375 // + faces sorted to correspond
2376 // e.g.
2377 // owner | neighbour
2378 // 0 | 2
2379 // 0 | 23
2380 // 0 | 71
2381 // 1 | 23
2382 // 1 | 24
2383 // 1 | 91
2384
2385 List<labelPair> ownerNeighbourPair(owner.size());
2386
2387 forAll(ownerNeighbourPair, oNI)
2388 {
2389 ownerNeighbourPair[oNI] = labelPair(owner[oNI], neighbour[oNI]);
2390 }
2391
2392 Info<< nl
2393 << "Sorting faces, owner and neighbour into upper triangular order"
2394 << endl;
2395
2396 labelList oldToNew(sortedOrder(ownerNeighbourPair));
2397 oldToNew = invert(oldToNew.size(), oldToNew);
2398
2399 inplaceReorder(oldToNew, faces);
2400 inplaceReorder(oldToNew, owner);
2401 inplaceReorder(oldToNew, neighbour);
2402}
2403
2404
2405void Foam::conformalVoronoiMesh::sortProcPatches
2406(
2407 List<DynamicList<face>>& patchFaces,
2408 List<DynamicList<label>>& patchOwners,
2409 List<DynamicList<label>>& patchPointPairSlaves,
2410 labelPairPairDynListList& patchSortingIndices
2411) const
2412{
2413 if (!Pstream::parRun())
2414 {
2415 return;
2416 }
2417
2418 forAll(patchSortingIndices, patchi)
2419 {
2420 faceList& faces = patchFaces[patchi];
2421 labelList& owner = patchOwners[patchi];
2422 DynamicList<label>& slaves = patchPointPairSlaves[patchi];
2423 DynamicList<labelPairPair>& sortingIndices
2424 = patchSortingIndices[patchi];
2425
2426 if (!sortingIndices.empty())
2427 {
2428 if
2429 (
2430 faces.size() != sortingIndices.size()
2431 || owner.size() != sortingIndices.size()
2432 || slaves.size() != sortingIndices.size()
2433 )
2434 {
2436 << "patch size and size of sorting indices is inconsistent "
2437 << " for patch " << patchi << nl
2438 << " faces.size() " << faces.size() << nl
2439 << " owner.size() " << owner.size() << nl
2440 << " slaves.size() " << slaves.size() << nl
2441 << " sortingIndices.size() "
2442 << sortingIndices.size()
2443 << exit(FatalError) << endl;
2444 }
2445
2446 labelList oldToNew(sortedOrder(sortingIndices));
2447 oldToNew = invert(oldToNew.size(), oldToNew);
2448
2449 inplaceReorder(oldToNew, sortingIndices);
2450 inplaceReorder(oldToNew, faces);
2451 inplaceReorder(oldToNew, owner);
2452 inplaceReorder(oldToNew, slaves);
2453 }
2454 }
2455}
2456
2457
2458void Foam::conformalVoronoiMesh::addPatches
2459(
2460 const label nInternalFaces,
2461 faceList& faces,
2462 labelList& owner,
2464 bitSet& boundaryFacesToRemove,
2465 const List<DynamicList<face>>& patchFaces,
2466 const List<DynamicList<label>>& patchOwners,
2467 const List<DynamicList<bool>>& indirectPatchFace
2468) const
2469{
2470 label nBoundaryFaces = 0;
2471
2472 forAll(patchFaces, p)
2473 {
2474 patchDicts[p].set("nFaces", patchFaces[p].size());
2475 patchDicts[p].set("startFace", nInternalFaces + nBoundaryFaces);
2476
2477 nBoundaryFaces += patchFaces[p].size();
2478 }
2479
2480 faces.setSize(nInternalFaces + nBoundaryFaces);
2481 owner.setSize(nInternalFaces + nBoundaryFaces);
2482 boundaryFacesToRemove.setSize(nInternalFaces + nBoundaryFaces);
2483
2484 label facei = nInternalFaces;
2485
2486 forAll(patchFaces, p)
2487 {
2488 forAll(patchFaces[p], f)
2489 {
2490 faces[facei] = patchFaces[p][f];
2491 owner[facei] = patchOwners[p][f];
2492 boundaryFacesToRemove[facei] = indirectPatchFace[p][f];
2493
2494 facei++;
2495 }
2496 }
2497}
2498
2499
2500void Foam::conformalVoronoiMesh::removeUnusedPoints
2501(
2502 faceList& faces,
2503 pointField& pts,
2504 labelList& boundaryPts
2505) const
2506{
2507 Info<< nl << "Removing unused points" << endl;
2508
2509 bitSet ptUsed(pts.size(), false);
2510
2511 // Scan all faces to find all of the points that are used
2512
2513 forAll(faces, fI)
2514 {
2515 const face& f = faces[fI];
2516
2517 ptUsed.set(f);
2518 }
2519
2520 label pointi = 0;
2521
2522 labelList oldToNew(pts.size(), label(-1));
2523
2524 // Move all of the used points to the start of the pointField and
2525 // truncate it
2526
2527 forAll(ptUsed, ptUI)
2528 {
2529 if (ptUsed.test(ptUI))
2530 {
2531 oldToNew[ptUI] = pointi++;
2532 }
2533 }
2534
2535 inplaceReorder(oldToNew, pts);
2536 inplaceReorder(oldToNew, boundaryPts);
2537
2538 Info<< " Removing "
2539 << returnReduce(pts.size() - pointi, sumOp<label>())
2540 << " unused points"
2541 << endl;
2542
2543 pts.setSize(pointi);
2544 boundaryPts.setSize(pointi);
2545
2546 // Renumber the faces to use the new point numbers
2547
2548 forAll(faces, fI)
2549 {
2550 inplaceRenumber(oldToNew, faces[fI]);
2551 }
2552}
2553
2554
2555Foam::labelList Foam::conformalVoronoiMesh::removeUnusedCells
2556(
2557 labelList& owner,
2558 labelList& neighbour
2559) const
2560{
2561 Info<< nl << "Removing unused cells" << endl;
2562
2563 bitSet cellUsed(vertexCount(), false);
2564
2565 // Scan all faces to find all of the cells that are used
2566
2567 cellUsed.set(owner);
2568 cellUsed.set(neighbour);
2569
2570 label celli = 0;
2571
2572 labelList oldToNew(cellUsed.size(), label(-1));
2573
2574 // Move all of the used cellCentres to the start of the pointField and
2575 // truncate it
2576
2577 forAll(cellUsed, cellUI)
2578 {
2579 if (cellUsed.test(cellUI))
2580 {
2581 oldToNew[cellUI] = celli++;
2582 }
2583 }
2584
2585 labelList newToOld(invert(celli, oldToNew));
2586
2587 // Find all of the unused cells, create a list of them, then
2588 // subtract one from each owner and neighbour entry for each of
2589 // the unused cell indices that it is above.
2590
2591 DynamicList<label> unusedCells;
2592
2593 forAll(cellUsed, cUI)
2594 {
2595 if (!cellUsed.test(cUI))
2596 {
2597 unusedCells.append(cUI);
2598 }
2599 }
2600
2601 if (unusedCells.size() > 0)
2602 {
2603 Info<< " Removing "
2604 << returnReduce(unusedCells.size(), sumOp<label>())
2605 << " unused cell labels" << endl;
2606
2607 forAll(owner, oI)
2608 {
2609 label& o = owner[oI];
2610
2611 o -= findLower(unusedCells, o) + 1;
2612 }
2613
2614 forAll(neighbour, nI)
2615 {
2616 label& n = neighbour[nI];
2617
2618 n -= findLower(unusedCells, n) + 1;
2619 }
2620 }
2621
2622 return newToOld;
2623}
2624
2625
2626// ************************************************************************* //
Various functions to operate on Lists.
label n
labelList faceLabels(nFaceLabels)
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
bool set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition List.H:469
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
void setSize(label n)
Alias for resize().
Definition List.H:536
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
void setSize(const label n)
Same as resize().
Definition PtrList.H:357
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
static const Form max
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
A collection of cell labels.
Definition cellSet.H:50
static void timeCheck(const Time &runTime, const string &description=string::null, const bool check=true)
Write the elapsedCpuTime and memory usage, with an optional.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
@ REGEX_RECURSIVE
Definition keyType.H:88
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces, const bool dryRun=false)
Check mesh with mesh settings in dict. Collects incorrect faces.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
static bool checkFaceDotProduct(const bool report, const scalar orthWarn, const polyMesh &, const vectorField &cellCentres, const vectorField &faceAreas, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
See primitiveMesh.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return pointer to a new patch created on freestore from components.
label nCells() const noexcept
Number of mesh cells.
Neighbour processor patch.
static word newName(const label myProcNo, const label neighbProcNo)
Return the name of a processorPolyPatch ("procBoundary..") constructed from the pair of processor IDs...
A class for handling words, derived from Foam::string.
Definition word.H:66
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
volScalarField & p
const label nProcPatches
const polyBoundaryMesh & patches
bool coupled
dynamicFvMesh & mesh
Foam::autoPtr< Foam::dynamicFvMesh > meshPtr
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
const pointField & points
const cellShapeList & cells
#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
void drawDelaunayCell(Ostream &os, const CellHandle &c, label offset=0)
Draws a tet cell to an output stream. The offset is supplied as the tet.
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
void writeProcessorInterface(const fileName &fName, const Triangulation &t, const faceList &faces)
Write the processor interface to an OBJ file.
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
Definition debug.C:45
const wordList surface
Standard surface field types (scalar, vector, tensor, etc).
const wordList internal
Standard volume internal field types (scalar, vector, tensor, etc).
const std::string patch
OpenFOAM patch number as a std::string.
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
List< word > wordList
List of word.
Definition fileName.H:60
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
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)
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition UListI.H:539
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
static const Identity< scalar > I
Definition Identity.H:100
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Binary search to find the index of the last element in a sorted list that is less than value.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition ListOps.C:28
Pair< labelPair > labelPairPair
Pair of labelPair.
Definition labelPair.H:32
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
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
Foam::point pointFromPoint
wordList patchNames(nPatches)
PtrList< dictionary > patchDicts
labelList f(nPoints)
label nPatches
dictionary dict
const pointField & pts
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299