Loading...
Searching...
No Matches
polyDualMesh.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2022,2025 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
27InClass
28 polyDualMesh
29
30\*---------------------------------------------------------------------------*/
31
32#include "polyDualMesh.H"
33#include "meshTools.H"
34#include "OFstream.H"
35#include "Time.H"
36#include "SortableList.H"
37#include "pointSet.H"
38#include "bitSet.H"
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42namespace Foam
43{
45}
46
47// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
49// Determine order for faces:
50// - upper-triangular order for internal faces
51// - external faces after internal faces (were already so)
52Foam::labelList Foam::polyDualMesh::getFaceOrder
53(
54 const labelList& faceOwner,
55 const labelList& faceNeighbour,
56 const cellList& cells,
57 label& nInternalFaces
58)
59{
60 labelList oldToNew(faceOwner.size(), -1);
61
62 // First unassigned face
63 label newFacei = 0;
64
65 forAll(cells, celli)
66 {
67 const labelList& cFaces = cells[celli];
68
69 SortableList<label> nbr(cFaces.size());
70
71 forAll(cFaces, i)
72 {
73 label facei = cFaces[i];
74
75 label nbrCelli = faceNeighbour[facei];
76
77 if (nbrCelli != -1)
78 {
79 // Internal face. Get cell on other side.
80 if (nbrCelli == celli)
81 {
82 nbrCelli = faceOwner[facei];
83 }
84
85 if (celli < nbrCelli)
86 {
87 // Celli is master
88 nbr[i] = nbrCelli;
89 }
90 else
91 {
92 // nbrCell is master. Let it handle this face.
93 nbr[i] = -1;
94 }
95 }
96 else
97 {
98 // External face. Do later.
99 nbr[i] = -1;
100 }
101 }
102
103 nbr.sort();
104
105 forAll(nbr, i)
106 {
107 if (nbr[i] != -1)
108 {
109 oldToNew[cFaces[nbr.indices()[i]]] = newFacei++;
110 }
111 }
112 }
113
114 nInternalFaces = newFacei;
115
116 Pout<< "nInternalFaces:" << nInternalFaces << endl;
117 Pout<< "nFaces:" << faceOwner.size() << endl;
118 Pout<< "nCells:" << cells.size() << endl;
119
120 // Leave patch faces intact.
121 for (label facei = newFacei; facei < faceOwner.size(); facei++)
122 {
123 oldToNew[facei] = facei;
124 }
125
126
127 // Check done all faces.
128 forAll(oldToNew, facei)
129 {
130 if (oldToNew[facei] == -1)
131 {
133 << "Did not determine new position"
134 << " for face " << facei
135 << abort(FatalError);
136 }
137 }
138
139 return oldToNew;
140}
141
142
143// Get the two edges on facei using pointi. Returns them such that the order
144// - otherVertex of e0
145// - pointi
146// - otherVertex(pointi) of e1
147// is in face order
148void Foam::polyDualMesh::getPointEdges
149(
150 const primitivePatch& patch,
151 const label facei,
152 const label pointi,
153 label& e0,
154 label& e1
155)
156{
157 const labelList& fEdges = patch.faceEdges()[facei];
158 const face& f = patch.localFaces()[facei];
159
160 e0 = -1;
161 e1 = -1;
162
163 forAll(fEdges, i)
164 {
165 label edgeI = fEdges[i];
166
167 const edge& e = patch.edges()[edgeI];
168
169 if (e[0] == pointi)
170 {
171 // One of the edges using pointi. Check which one.
172 label index = f.find(pointi);
173
174 if (f.nextLabel(index) == e[1])
175 {
176 e1 = edgeI;
177 }
178 else
179 {
180 e0 = edgeI;
181 }
182
183 if (e0 != -1 && e1 != -1)
184 {
185 return;
186 }
187 }
188 else if (e[1] == pointi)
189 {
190 // One of the edges using pointi. Check which one.
191 label index = f.find(pointi);
192
193 if (f.nextLabel(index) == e[0])
194 {
195 e1 = edgeI;
196 }
197 else
198 {
199 e0 = edgeI;
200 }
201
202 if (e0 != -1 && e1 != -1)
203 {
204 return;
205 }
206 }
207 }
208
210 << " vertices:" << patch.localFaces()[facei]
211 << " that uses point:" << pointi
212 << abort(FatalError);
213}
214
215
216// Collect the face on an external point of the patch.
217Foam::labelList Foam::polyDualMesh::collectPatchSideFace
218(
219 const polyPatch& patch,
220 const label patchToDualOffset,
221 const labelList& edgeToDualPoint,
222 const labelList& pointToDualPoint,
223 const label pointi,
224
225 label& edgeI
226)
227{
228 // Construct face by walking around e[eI] starting from
229 // patchEdgeI
230
231 label meshPointi = patch.meshPoints()[pointi];
232 const labelList& pFaces = patch.pointFaces()[pointi];
233
234 DynamicList<label> dualFace;
235
236 if (pointToDualPoint[meshPointi] >= 0)
237 {
238 // Number of pFaces + 2 boundary edge + feature point
239 dualFace.setCapacity(pFaces.size()+2+1);
240 // Store dualVertex for feature edge
241 dualFace.append(pointToDualPoint[meshPointi]);
242 }
243 else
244 {
245 dualFace.setCapacity(pFaces.size()+2);
246 }
247
248 // Store dual vertex for starting edge.
249 if (edgeToDualPoint[patch.meshEdges()[edgeI]] < 0)
250 {
252 << abort(FatalError);
253 }
254
255 dualFace.append(edgeToDualPoint[patch.meshEdges()[edgeI]]);
256
257 label facei = patch.edgeFaces()[edgeI][0];
258
259 // Check order of vertices of edgeI in face to see if we need to reverse.
260 bool reverseFace;
261
262 label e0, e1;
263 getPointEdges(patch, facei, pointi, e0, e1);
264
265 if (e0 == edgeI)
266 {
267 reverseFace = true;
268 }
269 else
270 {
271 reverseFace = false;
272 }
273
274 while (true)
275 {
276 // Store dual vertex for facei.
277 dualFace.append(facei + patchToDualOffset);
278
279 // Cross face to other edge on pointi
280 label e0, e1;
281 getPointEdges(patch, facei, pointi, e0, e1);
282
283 if (e0 == edgeI)
284 {
285 edgeI = e1;
286 }
287 else
288 {
289 edgeI = e0;
290 }
291
292 if (edgeToDualPoint[patch.meshEdges()[edgeI]] >= 0)
293 {
294 // Feature edge. Insert dual vertex for edge.
295 dualFace.append(edgeToDualPoint[patch.meshEdges()[edgeI]]);
296 }
297
298 const labelList& eFaces = patch.edgeFaces()[edgeI];
299
300 if (eFaces.size() == 1)
301 {
302 // Reached other edge of patch
303 break;
304 }
305
306 // Cross edge to other face.
307 if (eFaces[0] == facei)
308 {
309 facei = eFaces[1];
310 }
311 else
312 {
313 facei = eFaces[0];
314 }
315 }
316
317 dualFace.shrink();
318
319 if (reverseFace)
320 {
321 reverse(dualFace);
322 }
323
324 return dualFace;
325}
326
327
328// Collect face around pointi which is not on the outside of the patch.
329// Returns the vertices of the face and the indices in these vertices of
330// any points which are on feature edges.
331void Foam::polyDualMesh::collectPatchInternalFace
332(
333 const polyPatch& patch,
334 const label patchToDualOffset,
335 const labelList& edgeToDualPoint,
336
337 const label pointi,
338 const label startEdgeI,
339
340 labelList& dualFace2,
341 labelList& featEdgeIndices2
342)
343{
344 // Construct face by walking around pointi starting from startEdgeI
345 const labelList& meshEdges = patch.meshEdges();
346 const labelList& pFaces = patch.pointFaces()[pointi];
347
348 // Vertices of dualFace
349 DynamicList<label> dualFace(pFaces.size());
350 // Indices in dualFace of vertices added because of feature edge
351 DynamicList<label> featEdgeIndices(dualFace.size());
352
353
354 label edgeI = startEdgeI;
355 label facei = patch.edgeFaces()[edgeI][0];
356
357 // Check order of vertices of edgeI in face to see if we need to reverse.
358 bool reverseFace;
359
360 label e0, e1;
361 getPointEdges(patch, facei, pointi, e0, e1);
362
363 if (e0 == edgeI)
364 {
365 reverseFace = true;
366 }
367 else
368 {
369 reverseFace = false;
370 }
371
372 while (true)
373 {
374 // Insert dual vertex for face
375 dualFace.append(facei + patchToDualOffset);
376
377 // Cross face to other edge on pointi
378 label e0, e1;
379 getPointEdges(patch, facei, pointi, e0, e1);
380
381 if (e0 == edgeI)
382 {
383 edgeI = e1;
384 }
385 else
386 {
387 edgeI = e0;
388 }
389
390 if (edgeToDualPoint[meshEdges[edgeI]] >= 0)
391 {
392 // Feature edge. Insert dual vertex for edge.
393 dualFace.append(edgeToDualPoint[meshEdges[edgeI]]);
394 featEdgeIndices.append(dualFace.size()-1);
395 }
396
397 if (edgeI == startEdgeI)
398 {
399 break;
400 }
401
402 // Cross edge to other face.
403 const labelList& eFaces = patch.edgeFaces()[edgeI];
404
405 if (eFaces[0] == facei)
406 {
407 facei = eFaces[1];
408 }
409 else
410 {
411 facei = eFaces[0];
412 }
413 }
414
415 dualFace2.transfer(dualFace);
416
417 featEdgeIndices2.transfer(featEdgeIndices);
418
419 if (reverseFace)
420 {
421 reverse(dualFace2);
422
423 // Correct featEdgeIndices for change in dualFace2
424 forAll(featEdgeIndices2, i)
425 {
426 featEdgeIndices2[i] = dualFace2.size() -1 - featEdgeIndices2[i];
427 }
428 // Reverse indices (might not be necessary but do anyway)
429 reverse(featEdgeIndices2);
430 }
431}
432
433
434void Foam::polyDualMesh::splitFace
435(
436 const polyPatch& patch,
437 const labelList& pointToDualPoint,
438
439 const label pointi,
440 const labelList& dualFace,
441 const labelList& featEdgeIndices,
442
443 DynamicList<face>& dualFaces,
444 DynamicList<label>& dualOwner,
445 DynamicList<label>& dualNeighbour,
446 DynamicList<label>& dualRegion
447)
448{
449
450 // Split face because of feature edges/points
451 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
452 label meshPointi = patch.meshPoints()[pointi];
453
454 if (pointToDualPoint[meshPointi] >= 0)
455 {
456 // Feature point. Do face-centre decomposition.
457
458 if (featEdgeIndices.size() < 2)
459 {
460 // Feature point but no feature edges. Not handled at the moment
461 dualFaces.append(face(dualFace));
462 dualOwner.append(meshPointi);
463 dualNeighbour.append(-1);
464 dualRegion.append(patch.index());
465 }
466 else
467 {
468 // Do 'face-centre' decomposition. Start from first feature
469 // edge create face up until next feature edge.
470
471 forAll(featEdgeIndices, i)
472 {
473 label startFp = featEdgeIndices[i];
474
475 label endFp = featEdgeIndices[(i+1) % featEdgeIndices.size()];
476
477 // Collect face from startFp to endFp
478 label sz = 0;
479
480 if (endFp > startFp)
481 {
482 sz = endFp - startFp + 2;
483 }
484 else
485 {
486 sz = endFp + dualFace.size() - startFp + 2;
487 }
488 face subFace(sz);
489
490 // feature point becomes face centre.
491 subFace[0] = pointToDualPoint[patch.meshPoints()[pointi]];
492
493 // Copy from startFp up to endFp.
494 for (label subFp = 1; subFp < subFace.size(); subFp++)
495 {
496 subFace[subFp] = dualFace[startFp];
497
498 startFp = (startFp+1) % dualFace.size();
499 }
500
501 dualFaces.append(face(subFace));
502 dualOwner.append(meshPointi);
503 dualNeighbour.append(-1);
504 dualRegion.append(patch.index());
505 }
506 }
507 }
508 else
509 {
510 // No feature point. Check if feature edges.
511 if (featEdgeIndices.size() < 2)
512 {
513 // Not enough feature edges. No split.
514 dualFaces.append(face(dualFace));
515 dualOwner.append(meshPointi);
516 dualNeighbour.append(-1);
517 dualRegion.append(patch.index());
518 }
519 else
520 {
521 // Multiple feature edges but no feature point.
522 // Split face along feature edges. Gives weird result if
523 // number of feature edges > 2.
524
525 // Storage for new face
526 DynamicList<label> subFace(dualFace.size());
527
528 forAll(featEdgeIndices, featI)
529 {
530 label startFp = featEdgeIndices[featI];
531 label endFp = featEdgeIndices[featEdgeIndices.fcIndex(featI)];
532
533 label fp = startFp;
534
535 while (true)
536 {
537 label vertI = dualFace[fp];
538
539 subFace.append(vertI);
540
541 if (fp == endFp)
542 {
543 break;
544 }
545
546 fp = dualFace.fcIndex(fp);
547 }
548
549 if (subFace.size() > 2)
550 {
551 // Enough vertices to create a face from.
552 dualFaces.append(face(subFace));
553 dualOwner.append(meshPointi);
554 dualNeighbour.append(-1);
555 dualRegion.append(patch.index());
556
557 subFace.clear();
558 }
559 }
560 // Check final face.
561 if (subFace.size() > 2)
562 {
563 // Enough vertices to create a face from.
564 dualFaces.append(face(subFace));
565 dualOwner.append(meshPointi);
566 dualNeighbour.append(-1);
567 dualRegion.append(patch.index());
568
569 subFace.clear();
570 }
571 }
572 }
573}
574
575
576// Create boundary face for every point in patch
577void Foam::polyDualMesh::dualPatch
578(
579 const polyPatch& patch,
580 const label patchToDualOffset,
581 const labelList& edgeToDualPoint,
582 const labelList& pointToDualPoint,
583
584 const pointField& dualPoints,
585
586 DynamicList<face>& dualFaces,
587 DynamicList<label>& dualOwner,
588 DynamicList<label>& dualNeighbour,
589 DynamicList<label>& dualRegion
590)
591{
592 const labelList& meshEdges = patch.meshEdges();
593
594 // Whether edge has been done.
595 // 0 : not
596 // 1 : done e.start()
597 // 2 : done e.end()
598 // 3 : done both
599 labelList doneEdgeSide(meshEdges.size(), Zero);
600
601 bitSet donePoint(patch.nPoints(), false);
602
603
604 // Do points on edge of patch
605 // ~~~~~~~~~~~~~~~~~~~~~~~~~~
606
607 forAll(doneEdgeSide, patchEdgeI)
608 {
609 const labelList& eFaces = patch.edgeFaces()[patchEdgeI];
610
611 if (eFaces.size() == 1)
612 {
613 const edge& e = patch.edges()[patchEdgeI];
614
615 forAll(e, eI)
616 {
617 label bitMask = 1<<eI;
618
619 if ((doneEdgeSide[patchEdgeI] & bitMask) == 0)
620 {
621 // Construct face by walking around e[eI] starting from
622 // patchEdgeI
623
624 label pointi = e[eI];
625
626 label edgeI = patchEdgeI;
627 labelList dualFace
628 (
629 collectPatchSideFace
630 (
631 patch,
632 patchToDualOffset,
633 edgeToDualPoint,
634 pointToDualPoint,
635
636 pointi,
637 edgeI
638 )
639 );
640
641 // Now edgeI is end edge. Mark as visited
642 if (patch.edges()[edgeI][0] == pointi)
643 {
644 doneEdgeSide[edgeI] |= 1;
645 }
646 else
647 {
648 doneEdgeSide[edgeI] |= 2;
649 }
650
651 dualFaces.append(face(dualFace));
652 dualOwner.append(patch.meshPoints()[pointi]);
653 dualNeighbour.append(-1);
654 dualRegion.append(patch.index());
655
656 doneEdgeSide[patchEdgeI] |= bitMask;
657 donePoint.set(pointi);
658 }
659 }
660 }
661 }
662
663
664
665 // Do patch-internal points
666 // ~~~~~~~~~~~~~~~~~~~~~~~~
667
668 forAll(donePoint, pointi)
669 {
670 if (!donePoint.test(pointi))
671 {
672 labelList dualFace, featEdgeIndices;
673
674 collectPatchInternalFace
675 (
676 patch,
677 patchToDualOffset,
678 edgeToDualPoint,
679 pointi,
680 patch.pointEdges()[pointi][0], // Arbitrary starting edge
681
682 dualFace,
683 featEdgeIndices
684 );
685
686 //- Either keep in one piece or split face according to feature.
687
689 //dualFaces.append(face(dualFace));
690 //dualOwner.append(patch.meshPoints()[pointi]);
691 //dualNeighbour.append(-1);
692 //dualRegion.append(patch.index());
693
694 splitFace
695 (
696 patch,
697 pointToDualPoint,
698 pointi,
699 dualFace,
700 featEdgeIndices,
701
702 dualFaces,
703 dualOwner,
704 dualNeighbour,
705 dualRegion
706 );
707
708 donePoint[pointi] = true;
709 }
710 }
711}
712
713
714void Foam::polyDualMesh::calcDual
715(
716 const polyMesh& mesh,
717 const labelList& featureEdges,
718 const labelList& featurePoints
719)
720{
722 const label nIntFaces = mesh.nInternalFaces();
723
724
725 // Get patch for all of outside
726 primitivePatch allBoundary(patches.faces(), mesh.points());
727
728 // Correspondence from patch edge to mesh edge.
729 labelList meshEdges
730 (
731 allBoundary.meshEdges
732 (
733 mesh.edges(),
735 )
736 );
737
738 {
739 pointSet nonManifoldPoints
740 (
741 mesh,
742 "nonManifoldPoints",
743 meshEdges.size() / 100
744 );
745
746 allBoundary.checkPointManifold(true, &nonManifoldPoints);
747
748 if (nonManifoldPoints.size())
749 {
750 nonManifoldPoints.write();
751
753 << "There are " << nonManifoldPoints.size() << " points where"
754 << " the outside of the mesh is non-manifold." << nl
755 << "Such a mesh cannot be converted to a dual." << nl
756 << "Writing points at non-manifold sites to pointSet "
757 << nonManifoldPoints.name()
758 << exit(FatalError);
759 }
760 }
761
762
763 // Assign points
764 // ~~~~~~~~~~~~~
765
766 // mesh label dualMesh vertex
767 // ---------- ---------------
768 // celli celli
769 // facei nCells+facei-nIntFaces
770 // featEdgeI nCells+nFaces-nIntFaces+featEdgeI
771 // featPointi nCells+nFaces-nIntFaces+nFeatEdges+featPointi
772
773 pointField dualPoints
774 (
775 mesh.nCells() // cell centres
776 + mesh.nBoundaryFaces() // boundary face centres
777 + featureEdges.size() // additional boundary edges
778 + featurePoints.size() // additional boundary points
779 );
780
781 label dualPointi = 0;
782
783
784 // Cell centres.
785 const pointField& cellCentres = mesh.cellCentres();
786
787 cellPoint_.setSize(cellCentres.size());
788
789 forAll(cellCentres, celli)
790 {
791 cellPoint_[celli] = dualPointi;
792 dualPoints[dualPointi++] = cellCentres[celli];
793 }
794
795 // Boundary faces centres
796 const pointField& faceCentres = mesh.faceCentres();
797
798 boundaryFacePoint_.setSize(mesh.nBoundaryFaces());
799
800 for (label facei = nIntFaces; facei < mesh.nFaces(); facei++)
801 {
802 boundaryFacePoint_[facei - nIntFaces] = dualPointi;
803 dualPoints[dualPointi++] = faceCentres[facei];
804 }
805
806 // Edge status:
807 // >0 : dual point label (edge is feature edge)
808 // -1 : is boundary edge.
809 // -2 : is internal edge.
810 labelList edgeToDualPoint(mesh.nEdges(), -2);
811
812 forAll(meshEdges, patchEdgeI)
813 {
814 label edgeI = meshEdges[patchEdgeI];
815
816 edgeToDualPoint[edgeI] = -1;
817 }
818
819 forAll(featureEdges, i)
820 {
821 label edgeI = featureEdges[i];
822
823 const edge& e = mesh.edges()[edgeI];
824
825 edgeToDualPoint[edgeI] = dualPointi;
826 dualPoints[dualPointi++] = e.centre(mesh.points());
827 }
828
829
830
831 // Point status:
832 // >0 : dual point label (point is feature point)
833 // -1 : is point on edge of patch
834 // -2 : is point on patch (but not on edge)
835 // -3 : is internal point.
836 labelList pointToDualPoint(mesh.nPoints(), -3);
837
838 forAll(patches, patchi)
839 {
840 const labelList& meshPoints = patches[patchi].meshPoints();
841
842 forAll(meshPoints, i)
843 {
844 pointToDualPoint[meshPoints[i]] = -2;
845 }
846
847 const labelListList& loops = patches[patchi].edgeLoops();
848
849 forAll(loops, i)
850 {
851 const labelList& loop = loops[i];
852
853 forAll(loop, j)
854 {
855 pointToDualPoint[meshPoints[loop[j]]] = -1;
856 }
857 }
858 }
859
860 forAll(featurePoints, i)
861 {
862 label pointi = featurePoints[i];
863
864 pointToDualPoint[pointi] = dualPointi;
865 dualPoints[dualPointi++] = mesh.points()[pointi];
866 }
867
868
869 // Storage for new faces.
870 // Dynamic sized since we don't know size.
871
872 DynamicList<face> dynDualFaces(mesh.nEdges());
873 DynamicList<label> dynDualOwner(mesh.nEdges());
874 DynamicList<label> dynDualNeighbour(mesh.nEdges());
875 DynamicList<label> dynDualRegion(mesh.nEdges());
876
877
878 // Generate faces from edges on the boundary
879 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
880
881 forAll(meshEdges, patchEdgeI)
882 {
883 label edgeI = meshEdges[patchEdgeI];
884
885 const edge& e = mesh.edges()[edgeI];
886
887 label owner = -1;
888 label neighbour = -1;
889
890 if (e[0] < e[1])
891 {
892 owner = e[0];
893 neighbour = e[1];
894 }
895 else
896 {
897 owner = e[1];
898 neighbour = e[0];
899 }
900
901 // Find the boundary faces using the edge.
902 const labelList& patchFaces = allBoundary.edgeFaces()[patchEdgeI];
903
904 if (patchFaces.size() != 2)
905 {
907 << "Cannot handle edges with " << patchFaces.size()
908 << " connected boundary faces."
909 << abort(FatalError);
910 }
911
912 label face0 = patchFaces[0] + nIntFaces;
913 const face& f0 = mesh.faces()[face0];
914
915 label face1 = patchFaces[1] + nIntFaces;
916
917 // We want to start walking from patchFaces[0] or patchFaces[1],
918 // depending on which one uses owner,neighbour in the right order.
919
920 label startFacei = -1;
921 label endFacei = -1;
922
923 label index = f0.find(neighbour);
924
925 if (f0.nextLabel(index) == owner)
926 {
927 startFacei = face0;
928 endFacei = face1;
929 }
930 else
931 {
932 startFacei = face1;
933 endFacei = face0;
934 }
935
936 // Now walk from startFacei to cell to face to cell etc. until endFacei
937
938 DynamicList<label> dualFace;
939
940 if (edgeToDualPoint[edgeI] >= 0)
941 {
942 // Number of cells + 2 boundary faces + feature edge point
943 dualFace.setCapacity(mesh.edgeCells()[edgeI].size()+2+1);
944 // Store dualVertex for feature edge
945 dualFace.append(edgeToDualPoint[edgeI]);
946 }
947 else
948 {
949 dualFace.setCapacity(mesh.edgeCells()[edgeI].size()+2);
950 }
951
952 // Store dual vertex for starting face.
953 dualFace.append(mesh.nCells() + startFacei - nIntFaces);
954
955 label celli = mesh.faceOwner()[startFacei];
956 label facei = startFacei;
957
958 while (true)
959 {
960 // Store dual vertex from celli.
961 dualFace.append(celli);
962
963 // Cross cell to other face on edge.
964 label f0, f1;
965 meshTools::getEdgeFaces(mesh, celli, edgeI, f0, f1);
966
967 if (f0 == facei)
968 {
969 facei = f1;
970 }
971 else
972 {
973 facei = f0;
974 }
975
976 // Cross face to other cell.
977 if (facei == endFacei)
978 {
979 break;
980 }
981
982 if (mesh.faceOwner()[facei] == celli)
983 {
984 celli = mesh.faceNeighbour()[facei];
985 }
986 else
987 {
988 celli = mesh.faceOwner()[facei];
989 }
990 }
991
992 // Store dual vertex for endFace.
993 dualFace.append(mesh.nCells() + endFacei - nIntFaces);
994
995 dynDualFaces.append(face(dualFace.shrink()));
996 dynDualOwner.append(owner);
997 dynDualNeighbour.append(neighbour);
998 dynDualRegion.append(-1);
999
1000 {
1001 // Check orientation.
1002 const face& f = dynDualFaces.last();
1003 const vector areaNorm = f.areaNormal(dualPoints);
1004
1005 if (((mesh.points()[owner] - dualPoints[f[0]]) & areaNorm) > 0)
1006 {
1008 << " on boundary edge:" << edgeI
1009 << mesh.points()[mesh.edges()[edgeI][0]]
1010 << mesh.points()[mesh.edges()[edgeI][1]]
1011 << endl;
1012 }
1013 }
1014 }
1015
1016
1017 // Generate faces from internal edges
1018 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1019
1020 forAll(edgeToDualPoint, edgeI)
1021 {
1022 if (edgeToDualPoint[edgeI] == -2)
1023 {
1024 // Internal edge.
1025
1026 // Find dual owner, neighbour
1027
1028 const edge& e = mesh.edges()[edgeI];
1029
1030 label owner = -1;
1031 label neighbour = -1;
1032
1033 if (e[0] < e[1])
1034 {
1035 owner = e[0];
1036 neighbour = e[1];
1037 }
1038 else
1039 {
1040 owner = e[1];
1041 neighbour = e[0];
1042 }
1043
1044 // Get a starting cell
1045 const labelList& eCells = mesh.edgeCells()[edgeI];
1046
1047 label celli = eCells[0];
1048
1049 // Get the two faces on the cell and edge.
1050 label face0, face1;
1051 meshTools::getEdgeFaces(mesh, celli, edgeI, face0, face1);
1052
1053 // Find the starting face by looking at the order in which
1054 // the face uses the owner, neighbour
1055 const face& f0 = mesh.faces()[face0];
1056
1057 label index = f0.find(neighbour);
1058
1059 bool f0OrderOk = (f0.nextLabel(index) == owner);
1060
1061 label startFacei = -1;
1062
1063 if (f0OrderOk == (mesh.faceOwner()[face0] == celli))
1064 {
1065 startFacei = face0;
1066 }
1067 else
1068 {
1069 startFacei = face1;
1070 }
1071
1072 // Walk face-cell-face until starting face reached.
1073 DynamicList<label> dualFace(mesh.edgeCells()[edgeI].size());
1074
1075 label facei = startFacei;
1076
1077 while (true)
1078 {
1079 // Store dual vertex from celli.
1080 dualFace.append(celli);
1081
1082 // Cross cell to other face on edge.
1083 label f0, f1;
1084 meshTools::getEdgeFaces(mesh, celli, edgeI, f0, f1);
1085
1086 if (f0 == facei)
1087 {
1088 facei = f1;
1089 }
1090 else
1091 {
1092 facei = f0;
1093 }
1094
1095 // Cross face to other cell.
1096 if (facei == startFacei)
1097 {
1098 break;
1099 }
1100
1101 if (mesh.faceOwner()[facei] == celli)
1102 {
1103 celli = mesh.faceNeighbour()[facei];
1104 }
1105 else
1106 {
1107 celli = mesh.faceOwner()[facei];
1108 }
1109 }
1110
1111 dynDualFaces.append(face(dualFace.shrink()));
1112 dynDualOwner.append(owner);
1113 dynDualNeighbour.append(neighbour);
1114 dynDualRegion.append(-1);
1115
1116 {
1117 // Check orientation.
1118 const face& f = dynDualFaces.last();
1119 const vector areaNorm = f.areaNormal(dualPoints);
1120
1121 if (((mesh.points()[owner] - dualPoints[f[0]]) & areaNorm) > 0)
1122 {
1124 << " on internal edge:" << edgeI
1125 << mesh.points()[mesh.edges()[edgeI][0]]
1126 << mesh.points()[mesh.edges()[edgeI][1]]
1127 << endl;
1128 }
1129 }
1130 }
1131 }
1132
1133 // Dump faces.
1134 if (debug)
1135 {
1136 dynDualFaces.shrink();
1137 dynDualOwner.shrink();
1138 dynDualNeighbour.shrink();
1139 dynDualRegion.shrink();
1140
1141 OFstream str("dualInternalFaces.obj");
1142 Pout<< "polyDualMesh::calcDual : dumping internal faces to "
1143 << str.name() << endl;
1144
1145 forAll(dualPoints, dualPointi)
1146 {
1147 meshTools::writeOBJ(str, dualPoints[dualPointi]);
1148 }
1149 forAll(dynDualFaces, dualFacei)
1150 {
1151 const face& f = dynDualFaces[dualFacei];
1152
1153 str<< 'f';
1154 forAll(f, fp)
1155 {
1156 str<< ' ' << f[fp]+1;
1157 }
1158 str<< nl;
1159 }
1160 }
1161
1162 const label nInternalFaces = dynDualFaces.size();
1163
1164 // Outside faces
1165 // ~~~~~~~~~~~~~
1166
1167 forAll(patches, patchi)
1168 {
1169 const polyPatch& pp = patches[patchi];
1170
1171 dualPatch
1172 (
1173 pp,
1174 mesh.nCells() + pp.start() - nIntFaces,
1175 edgeToDualPoint,
1176 pointToDualPoint,
1177
1178 dualPoints,
1179
1180 dynDualFaces,
1181 dynDualOwner,
1182 dynDualNeighbour,
1183 dynDualRegion
1184 );
1185 }
1186
1187
1188 // Transfer face info to straight lists
1189 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1190 faceList dualFaces(std::move(dynDualFaces));
1191 labelList dualOwner(std::move(dynDualOwner));
1192 labelList dualNeighbour(std::move(dynDualNeighbour));
1193 labelList dualRegion(std::move(dynDualRegion));
1194
1195
1196 // Dump faces.
1197 if (debug)
1198 {
1199 OFstream str("dualFaces.obj");
1200 Pout<< "polyDualMesh::calcDual : dumping all faces to "
1201 << str.name() << endl;
1202
1203 forAll(dualPoints, dualPointi)
1204 {
1205 meshTools::writeOBJ(str, dualPoints[dualPointi]);
1206 }
1207 forAll(dualFaces, dualFacei)
1208 {
1209 const face& f = dualFaces[dualFacei];
1210
1211 str<< 'f';
1212 forAll(f, fp)
1213 {
1214 str<< ' ' << f[fp]+1;
1215 }
1216 str<< nl;
1217 }
1218 }
1219
1220 // Create cells.
1221 cellList dualCells(mesh.nPoints());
1222
1223 // unnecessary...
1224 // forAll(dualCells, celli)
1225 // {
1226 // dualCells[celli].clear();
1227 // }
1228
1229 forAll(dualOwner, facei)
1230 {
1231 label celli = dualOwner[facei];
1232
1233 if (celli != -1)
1234 {
1235 dualCells[celli].push_back(facei);
1236 }
1237 }
1238 forAll(dualNeighbour, facei)
1239 {
1240 label celli = dualNeighbour[facei];
1241
1242 if (celli != -1)
1243 {
1244 dualCells[celli].push_back(facei);
1245 }
1246 }
1247
1248
1249 // Do upper-triangular face ordering. Determines face reordering map and
1250 // number of internal faces.
1251 label dummy;
1252
1253 labelList oldToNew
1254 (
1255 getFaceOrder
1256 (
1257 dualOwner,
1258 dualNeighbour,
1259 dualCells,
1260 dummy
1261 )
1262 );
1263
1264 // Reorder faces.
1265 inplaceReorder(oldToNew, dualFaces);
1266 inplaceReorder(oldToNew, dualOwner);
1267 inplaceReorder(oldToNew, dualNeighbour);
1268 inplaceReorder(oldToNew, dualRegion);
1269 forAll(dualCells, celli)
1270 {
1271 inplaceRenumber(oldToNew, dualCells[celli]);
1272 }
1273
1274
1275 // Create patches
1276 labelList patchSizes(patches.size(), Zero);
1277
1278 forAll(dualRegion, facei)
1279 {
1280 if (dualRegion[facei] >= 0)
1281 {
1282 patchSizes[dualRegion[facei]]++;
1283 }
1284 }
1285
1286 labelList patchStarts(patches.size(), Zero);
1287
1288 label facei = nInternalFaces;
1289
1290 forAll(patches, patchi)
1291 {
1292 patchStarts[patchi] = facei;
1293
1294 facei += patchSizes[patchi];
1295 }
1296
1297
1298 Pout<< "nFaces:" << dualFaces.size()
1299 << " patchSizes:" << patchSizes
1300 << " patchStarts:" << patchStarts
1301 << endl;
1302
1303
1304 // Add patches. First add zero sized (since mesh still 0 size)
1305 polyPatchList dualPatches(patches.size());
1306
1307 forAll(patches, patchi)
1308 {
1309 const polyPatch& pp = patches[patchi];
1310
1311 dualPatches.set
1312 (
1313 patchi,
1314 pp.clone
1315 (
1316 boundaryMesh(),
1317 patchi,
1318 0, //patchSizes[patchi],
1319 0 //patchStarts[patchi]
1320 )
1321 );
1322 }
1323 addPatches(dualPatches);
1324
1325 // Assign to mesh.
1326 resetPrimitives
1327 (
1328 autoPtr<pointField>::New(std::move(dualPoints)),
1329 autoPtr<faceList>::New(std::move(dualFaces)),
1330 autoPtr<labelList>::New(std::move(dualOwner)),
1331 autoPtr<labelList>::New(std::move(dualNeighbour)),
1332 patchSizes,
1333 patchStarts
1334 );
1335}
1337
1338// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1339
1340// Construct from components
1341Foam::polyDualMesh::polyDualMesh(const IOobject& io)
1342:
1343 polyMesh(io),
1344 cellPoint_
1345 (
1346 IOobject
1347 (
1348 "cellPoint",
1349 time().findInstance(meshDir(), "cellPoint"),
1350 meshSubDir,
1351 *this,
1354 )
1355 ),
1356 boundaryFacePoint_
1357 (
1358 IOobject
1359 (
1360 "boundaryFacePoint",
1361 time().findInstance(meshDir(), "boundaryFacePoint"),
1362 meshSubDir,
1363 *this,
1366 )
1368{}
1369
1370
1371// Construct from polyMesh
1372Foam::polyDualMesh::polyDualMesh
1373(
1374 const polyMesh& mesh,
1375 const labelList& featureEdges,
1376 const labelList& featurePoints
1377)
1378:
1379 polyMesh(mesh, Zero),
1380 cellPoint_
1381 (
1382 IOobject
1383 (
1384 "cellPoint",
1385 time().findInstance(meshDir(), "faces"),
1386 meshSubDir,
1387 *this,
1388 IOobject::NO_READ,
1389 IOobject::AUTO_WRITE
1390 ),
1391 labelList(mesh.nCells(), -1)
1392 ),
1393 boundaryFacePoint_
1394 (
1395 IOobject
1396 (
1397 "boundaryFacePoint",
1398 time().findInstance(meshDir(), "faces"),
1399 meshSubDir,
1400 *this,
1401 IOobject::NO_READ,
1402 IOobject::AUTO_WRITE
1403 ),
1404 labelList(mesh.nBoundaryFaces(), -1)
1405 )
1406{
1407 calcDual(mesh, featureEdges, featurePoints);
1408}
1409
1410
1411// Construct from polyMesh and feature angle
1412Foam::polyDualMesh::polyDualMesh
1413(
1414 const polyMesh& mesh,
1415 const scalar featureCos
1416)
1417:
1418 polyMesh(mesh, Zero),
1419 cellPoint_
1420 (
1421 IOobject
1422 (
1423 "cellPoint",
1424 time().findInstance(meshDir(), "faces"),
1425 meshSubDir,
1426 *this,
1427 IOobject::NO_READ,
1428 IOobject::AUTO_WRITE
1429 ),
1430 labelList(mesh.nCells(), -1)
1431 ),
1432 boundaryFacePoint_
1433 (
1434 IOobject
1435 (
1436 "boundaryFacePoint",
1437 time().findInstance(meshDir(), "faces"),
1438 meshSubDir,
1439 *this,
1440 IOobject::NO_READ,
1441 IOobject::AUTO_WRITE
1442 ),
1443 labelList(mesh.nBoundaryFaces(), -1)
1444 )
1445{
1446 labelList featureEdges, featurePoints;
1447
1448 calcFeatures(mesh, featureCos, featureEdges, featurePoints);
1449 calcDual(mesh, featureEdges, featurePoints);
1450}
1451
1452
1454(
1455 const polyMesh& mesh,
1456 const scalar featureCos,
1457 labelList& featureEdges,
1458 labelList& featurePoints
1459)
1460{
1461 const polyBoundaryMesh& patches = mesh.boundaryMesh();
1462
1463 // Create big primitivePatch for all outside.
1464 primitivePatch allBoundary(patches.faces(), mesh.points());
1465
1466 // For ease of use store patch number per face in allBoundary.
1467 labelList allRegion(allBoundary.size());
1468
1469 forAll(patches, patchi)
1470 {
1471 const polyPatch& pp = patches[patchi];
1472
1473 allRegion.slice(pp.offset(), pp.size()) = patchi;
1474 }
1475
1476
1477 // Calculate patch/feature edges
1478 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1479
1480 const labelListList& edgeFaces = allBoundary.edgeFaces();
1481 const vectorField& faceNormals = allBoundary.faceNormals();
1482 const labelList& meshPoints = allBoundary.meshPoints();
1483
1484 bitSet isFeatureEdge(edgeFaces.size(), false);
1485
1486 forAll(edgeFaces, edgeI)
1487 {
1488 const labelList& eFaces = edgeFaces[edgeI];
1489
1490 if (eFaces.size() != 2)
1491 {
1492 // Non-manifold. Problem?
1493 const edge& e = allBoundary.edges()[edgeI];
1494
1496 << meshPoints[e[0]] << ' ' << meshPoints[e[1]]
1497 << " coords:" << mesh.points()[meshPoints[e[0]]]
1498 << mesh.points()[meshPoints[e[1]]]
1499 << " has more than 2 faces connected to it:"
1500 << eFaces.size() << endl;
1501
1502 isFeatureEdge.set(edgeI);
1503 }
1504 else if (allRegion[eFaces[0]] != allRegion[eFaces[1]])
1505 {
1506 isFeatureEdge.set(edgeI);
1507 }
1508 else if
1509 (
1510 (faceNormals[eFaces[0]] & faceNormals[eFaces[1]])
1511 < featureCos
1512 )
1513 {
1514 isFeatureEdge.set(edgeI);
1515 }
1516 }
1517
1518
1519 // Calculate feature points
1520 // ~~~~~~~~~~~~~~~~~~~~~~~~
1521
1522 const labelListList& pointEdges = allBoundary.pointEdges();
1523
1524 DynamicList<label> allFeaturePoints(pointEdges.size());
1525
1526 forAll(pointEdges, pointi)
1527 {
1528 const labelList& pEdges = pointEdges[pointi];
1529
1530 label nFeatEdges = 0;
1531
1532 forAll(pEdges, i)
1533 {
1534 if (isFeatureEdge.test(pEdges[i]))
1535 {
1536 ++nFeatEdges;
1537 }
1538 }
1539 if (nFeatEdges > 2)
1540 {
1541 // Store in mesh vertex label
1542 allFeaturePoints.append(allBoundary.meshPoints()[pointi]);
1543 }
1544 }
1545 featurePoints.transfer(allFeaturePoints);
1546
1547 if (debug)
1548 {
1549 OFstream str("featurePoints.obj");
1550
1551 Pout<< "polyDualMesh::calcFeatures : dumping feature points to "
1552 << str.name() << endl;
1553
1554 forAll(featurePoints, i)
1555 {
1556 label pointi = featurePoints[i];
1557 meshTools::writeOBJ(str, mesh.points()[pointi]);
1558 }
1559 }
1560
1561
1562 // Get all feature edges.
1563 labelList meshEdges
1564 (
1565 allBoundary.meshEdges
1566 (
1567 mesh.edges(),
1568 mesh.cellEdges(),
1570 (
1571 mesh.faceOwner(),
1572 allBoundary.size(),
1574 )
1575 )
1576 );
1577
1578 DynamicList<label> allFeatureEdges(isFeatureEdge.size());
1579 forAll(isFeatureEdge, edgeI)
1580 {
1581 if (isFeatureEdge.test(edgeI))
1582 {
1583 // Store in mesh edge label.
1584 allFeatureEdges.append(meshEdges[edgeI]);
1585 }
1586 }
1587 featureEdges.transfer(allFeatureEdges);
1589
1590
1591// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1592
1594{}
1595
1596
1597// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
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
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
label size() const noexcept
Number of entries.
Definition PackedList.H:392
const labelListList & pointEdges() const
Return point-edge addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using cell addressing.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
const labelListList & edgeFaces() const
Return edge-face addressing.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
SubList< T > slice(const label pos, label len=-1)
Return SubList slice (non-const access) - no range checking.
Definition SubList.H:258
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
T & last()
Access last element of the list, position [size()-1].
Definition UList.H:971
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
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
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition bitSet.H:334
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
const Time & time() const noexcept
Return time registry.
A set of point labels.
Definition pointSet.H:50
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const faceList::subList faces() const
Return mesh faces for the entire boundary.
Creates dual of polyMesh.
~polyDualMesh()
Destructor.
static void calcFeatures(const polyMesh &, const scalar featureCos, labelList &featureEdges, labelList &featurePoints)
Helper function to create feature edges and points based on.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir).
Definition polyMesh.C:826
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces).
const labelListList & pointEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const labelListList & cellEdges() const
const vectorField & faceCentres() const
label nInternalFaces() const noexcept
Number of internal faces.
const vectorField & cellCentres() const
label nPoints() const noexcept
Number of mesh points.
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
const labelListList & edgeFaces() const
const labelListList & edgeCells() const
const cellList & cells() const
label nEdges() const
Number of mesh edges.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
const cellShapeList & cells
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
void getEdgeFaces(const primitiveMesh &mesh, const label celli, const label edgeI, label &face0, label &face1)
Get faces on cell using edgeI. Throws error if no two found.
Definition meshTools.C:472
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition polyPatch.H:61
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
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition UListI.H:539
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
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
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
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...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};constexpr label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< DynamicList< face > > pFaces[nBCs]
labelList f(nPoints)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299