Loading...
Searching...
No Matches
tetDecomposer.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) 2015-2020,2022,2024 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
29#include "tetDecomposer.H"
30#include "meshTools.H"
31#include "polyMesh.H"
32#include "polyTopoChange.H"
33#include "mapPolyMesh.H"
34#include "OFstream.H"
35#include "edgeHashes.H"
36#include "syncTools.H"
38#include "triangle.H"
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42namespace Foam
43{
45}
46
47const Foam::Enum
48<
50>
52({
53 { decompositionType::FACE_CENTRE_TRIS, "faceCentre" },
54 { decompositionType::FACE_DIAG_TRIS, "faceDiagonal" },
55 { decompositionType::PYRAMID, "pyramid" },
56 { decompositionType::FACE_DIAG_QUADS, "faceDiagonalQuads" },
57});
58
59
60// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
61
62void Foam::tetDecomposer::modifyFace
63(
64 polyTopoChange& meshMod,
65 const face& f,
66 const label facei,
67 const label own,
68 const label nei,
69 const label patchi,
70 const label zoneI,
71 const bool zoneFlip
72) const
73{
74 if (own == nei)
75 {
77 << "Problem own:" << own
78 << " nei:" << nei << exit(FatalError);
79 }
80
81 // First usage of face. Modify.
82 if (nei == -1 || own < nei)
83 {
84 meshMod.modifyFace
85 (
86 f, // modified face
87 facei, // label of face
88 own, // owner
89 nei, // neighbour
90 false, // face flip
91 patchi, // patch for face
92 zoneI, // zone for face
93 zoneFlip // face flip in zone
94 );
95 }
96 else
97 {
98 meshMod.modifyFace
99 (
100 f.reverseFace(), // modified face
101 facei, // label of face
102 nei, // owner
103 own, // neighbour
104 true, // face flip
105 patchi, // patch for face
106 zoneI, // zone for face
107 !zoneFlip // face flip in zone
108 );
109 }
110}
111
112
113void Foam::tetDecomposer::addFace
114(
115 polyTopoChange& meshMod,
116 const face& f,
117 const label facei,
118 const label own,
119 const label nei,
120 const label masterPointID,
121 const label masterEdgeID,
122 const label masterFaceID,
123 const label patchi,
124 const label zoneI,
125 const bool zoneFlip
126) const
127{
128 if (own == nei)
129 {
131 << "Problem own:" << own
132 << " nei:" << nei << exit(FatalError);
133 }
134
135 // Second or more usage of face. Add.
136 if (nei == -1 || own < nei)
137 {
138 //const label newFacei =
139 meshMod.addFace
140 (
141 f, // modified face
142 own, // owner
143 nei, // neighbour
144 masterPointID, // master point
145 masterEdgeID, // master edge
146 masterFaceID, // master face
147 false, // face flip
148 patchi, // patch for face
149 zoneI, // zone for face
150 zoneFlip // face flip in zone
151 );
152 }
153 else
154 {
155 //const label newFacei =
156 meshMod.addFace
157 (
158 f.reverseFace(), // modified face
159 nei, // owner
160 own, // neighbour
161 masterPointID, // master point
162 masterEdgeID, // master edge
163 masterFaceID, // master face
164 true, // face flip
165 patchi, // patch for face
166 zoneI, // zone for face
167 !zoneFlip // face flip in zone
168 );
169 }
170}
171
172
173// Work out triangle index given the starting vertex in the face
174Foam::label Foam::tetDecomposer::triIndex(const label facei, const label fp)
175const
176{
177 const face& f = mesh_.faces()[facei];
178 const label fp0 = max(0, mesh_.tetBasePtIs()[facei]);
179
180 // Work out triangle index on this face. Assume 'fan' triangulation starting
181 // from fp0. E.g. 6 vertices on face -> 4 triangles. First and last triangle
182 // use consecutive vertices
183 //
184 // fp | triangle | vertices
185 // ------|----------|---------
186 // fp0 | 0 | fp0,fp0+1,fp0+2
187 // fp0+1 | 0 | ,,
188 // fp0+2 | 1 | fp0,fp0+2,fp0+3
189 // fp0+3 | 2 | fp0,fp0+3,fp0+4
190 // fp0+4 | 3 | fp0,fp0+4,fp0+3
191 // fp0+5 | 3 | ,,
192
193
194 // Work out triangle index on this face
195 label thisTrii;
196 if (fp == fp0)
197 {
198 thisTrii = 0;
199 }
200 else if (fp == f.fcIndex(fp0))
201 {
202 thisTrii = 0;
203 }
204 else if (fp == f.rcIndex(fp0))
205 {
206 thisTrii = f.size()-3;
207 }
208 else if (fp < fp0)
209 {
210 const label fpB = fp+f.size();
211 thisTrii = (fpB-fp0-1);
212 }
213 else
214 {
215 thisTrii = (fp-fp0-1);
216 }
217 return thisTrii;
218}
219
220
221void Foam::tetDecomposer::splitBoundaryFaces
222(
223 List<faceList>& boundaryQuads,
224 List<faceList>& boundaryTris
225) const
226{
227 // Work space
228 faceList quadFaces(1000);
229 faceList triFaces(1000);
230
231 const auto& pbm = mesh_.boundaryMesh();
232 for (const auto& pp : pbm)
233 {
234 if (pp.coupled() && refCast<const coupledPolyPatch>(pp).owner())
235 {
236 forAll(pp, i)
237 {
238 const label facei = pp.start()+i;
239 const face& meshFace = pp[i];
240
241 if (meshFace.size() > 4)
242 {
243 label trii = 0;
244 label quadi = 0;
245 meshFace.trianglesQuads
246 (
247 mesh_.points(),
248 trii,
249 quadi,
250 triFaces,
251 quadFaces
252 );
253
254 const label bFacei = facei-mesh_.nInternalFaces();
255
256 {
257 auto& faces = boundaryTris[bFacei];
258 faces.setSize(trii);
259 // Express as relative w.r.t. 0th vertex
260 for (label i = 0; i < trii; i++)
261 {
262 const auto& f = triFaces[i];
263 auto& verts = faces[i];
264 verts.setSize(f.size());
265 forAll(f, fp)
266 {
267 verts[fp] = meshFace.find(f[fp]);
268 }
269 }
270 }
271 {
272 auto& faces = boundaryQuads[bFacei];
273 faces.setSize(quadi);
274 // Express as relative w.r.t. 0th vertex
275 for (label i = 0; i < quadi; i++)
276 {
277 const auto& f = quadFaces[i];
278 auto& verts = faces[i];
279 verts.setSize(f.size());
280 forAll(f, fp)
281 {
282 verts[fp] = meshFace.find(f[fp]);
283 }
284 }
285 }
286 }
287 }
288 }
289 }
290
291 // Send coupled side indices to neighbour. Note: could also do re-indexing
292 // here...
294 (
295 mesh_,
296 boundaryTris,
297 [](faceList& result, const UList<face>& input)
298 {
299 if (result.empty())
300 {
301 result = input;
302 }
303 },
305 );
307 (
308 mesh_,
309 boundaryQuads,
310 [](faceList& result, const UList<face>& input)
311 {
312 if (result.empty())
313 {
314 result = input;
315 }
316 },
318 );
319}
320
321
322void Foam::tetDecomposer::relativeIndicesToFace
323(
324 const bool doFlip,
325 const face& meshFace,
326 const faceList& indexLists,
327 faceList& faces
328) const
329{
330 //faces.setSize(indexLists.size());
331
332 if (!doFlip)
333 {
334 forAll(indexLists, facei)
335 {
336 const auto& verts = indexLists[facei];
337 auto& f = faces[facei];
338 f.setSize(verts.size());
339
340 forAll(verts, fp)
341 {
342 f[fp] = meshFace[verts[fp]];
343 }
344 }
345 }
346 else
347 {
348 forAllReverse(indexLists, facei)
349 {
350 const auto& verts = indexLists[facei];
351 auto& f = faces[facei];
352 f.setSize(verts.size());
353
354 // - 0th vertex matches; walking order opposite
355 // - assemble in opposite order so as to flip normal
356 label destFp = verts.size()-1;
357 forAll(verts, fp)
358 {
359 if (verts[fp] == 0)
360 {
361 f[destFp] = meshFace[0];
362 }
363 else
364 {
365 f[destFp] = meshFace[meshFace.size()-verts[fp]];
366 }
367 --destFp;
368 }
369 }
370 }
371}
372
373
374void Foam::tetDecomposer::splitFace
375(
376 const List<faceList>& boundaryQuads,
377 const List<faceList>& boundaryTris,
378 const label facei,
379 const label patchi,
380 label& quadi,
381 faceList& quadFaces,
382 label& trii,
383 faceList& triFaces
384) const
385{
386 // Split face into quads (in quadFaces) and tris (in triFaces)
387
388 const face& f = mesh_.faces()[facei];
389
390 trii = 0;
391 quadi = 0;
392 if (patchi != -1 && mesh_.boundaryMesh()[patchi].coupled())
393 {
394 const auto& pp = mesh_.boundaryMesh()[patchi];
395 const bool owner =
397 const label bFacei = facei-mesh_.nInternalFaces();
398
399 // Triangles
400 trii = boundaryTris[bFacei].size();
401 relativeIndicesToFace
402 (
403 !owner,
404 f,
405 boundaryTris[bFacei],
406 triFaces
407 );
408
409 // Quads
410 quadi = boundaryQuads[bFacei].size();
411 relativeIndicesToFace
412 (
413 !owner,
414 f,
415 boundaryQuads[bFacei],
416 quadFaces
417 );
418 }
419 else if (f.size() == 4)
420 {
421 quadFaces[quadi++] = f;
422 }
423 else
424 {
425 f.trianglesQuads
426 (
427 mesh_.points(),
428 trii,
429 quadi,
430 triFaces,
431 quadFaces
432 );
433 }
434}
435
436
437void Foam::tetDecomposer::splitFace
438(
439 const List<faceList>& boundaryQuads,
440 const List<faceList>& boundaryTris,
441 const label facei,
442 const label patchi,
443 faceList& quadFaces,
444 faceList& triFaces,
445 faceList& subFaces
446) const
447{
448 const face& f = mesh_.faces()[facei];
449
450 if (f.size() <= 4)
451 {
452 subFaces.resize_nocopy(1);
453 subFaces[0] = f;
454 }
455 else
456 {
457 label trii;
458 label quadi;
459 splitFace
460 (
461 boundaryQuads,
462 boundaryTris,
463 facei,
464 patchi,
465 quadi,
466 quadFaces,
467 trii,
468 triFaces
469 );
470
471 // Collect into single face list
472 subFaces.resize_nocopy(quadi+trii);
473 {
474 label subFacei = 0;
475 for (label i = 0; i < quadi; i++)
476 {
477 subFaces[subFacei++] = quadFaces[i];
478 }
479 for (label i = 0; i < trii; i++)
480 {
481 subFaces[subFacei++] = triFaces[i];
482 }
484 }
485}
486
487
488// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
489
490Foam::tetDecomposer::tetDecomposer(const polyMesh& mesh)
492 mesh_(mesh)
493{}
494
495
496// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
497
499(
500 const decompositionType decomposeType,
501 const bitSet& decomposeCell,
502 polyTopoChange& meshMod
503)
504{
505 // Determine for every face whether it borders a cell that is decomposed
506 bitSet decomposeFace(mesh_.nFaces());
507 {
508 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
509 {
510 label own = mesh_.faceOwner()[facei];
511 label nei = mesh_.faceNeighbour()[facei];
512 if (decomposeCell[own] || decomposeCell[nei])
513 {
514 decomposeFace[facei] = true;
515 }
516 }
517
518 boolList neiDecomposeCell(mesh_.nBoundaryFaces());
519 forAll(neiDecomposeCell, bFacei)
520 {
521 label facei = mesh_.nInternalFaces()+bFacei;
522 label own = mesh_.faceOwner()[facei];
523 neiDecomposeCell[bFacei] = decomposeCell[own];
524 }
525 syncTools::swapBoundaryFaceList(mesh_, neiDecomposeCell);
526
527 for
528 (
529 label facei = mesh_.nInternalFaces();
530 facei < mesh_.nFaces();
531 facei++
532 )
533 {
534 label own = mesh_.faceOwner()[facei];
535 label bFacei = facei-mesh_.nInternalFaces();
536 if (decomposeCell[own] || neiDecomposeCell[bFacei])
537 {
538 decomposeFace[facei] = true;
540 }
541 }
542 setRefinement(decomposeType, decomposeCell, decomposeFace, meshMod);
543}
544
545
547(
548 const decompositionType decomposeType,
549 const bitSet& decomposeCell,
550 const bitSet& decomposeFace,
551 polyTopoChange& meshMod
552)
553{
554 if (debug)
555 {
556 // Check that current mesh makes sense
557 const auto& faces = mesh_.faces();
558 labelList nVerts(mesh_.nFaces());
559 forAll(faces, facei)
560 {
561 nVerts[facei] = faces[facei].size();
562 }
563 syncTools::swapFaceList(mesh_, nVerts);
564 forAll(nVerts, facei)
565 {
566 if (nVerts[facei] != faces[facei].size())
567 {
568 FatalErrorInFunction<< "problem with nVerts"
569 << exit(FatalError);
570 }
571 }
572 }
573 if (debug)
574 {
575 // Check that decomposeFace makes sense
576 bitSet newDecomposeFace(decomposeFace);
577 syncTools::swapFaceList(mesh_, newDecomposeFace);
578 forAll(newDecomposeFace, facei)
579 {
580 if (decomposeFace[facei] != newDecomposeFace[facei])
581 {
582 FatalErrorInFunction<< "problem with decomposeFace"
583 << exit(FatalError);
584 }
585 }
586 }
587
588
589 // For diagonal splitting : send over diagonals since two sides of
590 // coupled face might make different decisions
591 // TBD: probably also for FACE_DIAG_TRIS ?
592 faceList quadFaces;
593 faceList triFaces;
594 List<faceList> boundaryQuads;
595 List<faceList> boundaryTris;
596 faceList subFaces;
597
598 if (decomposeType == FACE_DIAG_QUADS)
599 {
600 // Pre-calculate coupled faces triangulation. Store as indices w.r.t.
601 // vertex 0.
602 // Note: could triangulate all faces here but trying to save some
603 // memory so only doing:
604 // - faces on proc boundaries
605 // - faces > 4 verts
606 quadFaces.resize_nocopy(1000);
607 triFaces.resize_nocopy(1000);
608 boundaryQuads.resize_nocopy(mesh_.nBoundaryFaces());
609 boundaryTris.resize_nocopy(mesh_.nBoundaryFaces());
610 splitBoundaryFaces(boundaryQuads, boundaryTris);
611 }
612
613
614 cellToPoint_.setSize(mesh_.nCells(), -1);
615 forAll(mesh_.cellCentres(), celli)
616 {
617 if (decomposeCell[celli])
618 {
619 // Any point on the cell
620 label masterPointi = mesh_.faces()[mesh_.cells()[celli][0]][0];
621
622 cellToPoint_[celli] = meshMod.addPoint
623 (
624 mesh_.cellCentres()[celli],
625 masterPointi,
626 -1,
627 true
628 );
629 }
630 }
631
632
633 // Add face centre points
634 if (decomposeType == FACE_CENTRE_TRIS)
635 {
636 faceToPoint_.setSize(mesh_.nFaces(), -1);
637 forAll(mesh_.faceCentres(), facei)
638 {
639 if (decomposeFace[facei])
640 {
641 // Any point on the face
642 const label masterPointi = mesh_.faces()[facei][0];
643
644 faceToPoint_[facei] = meshMod.addPoint
645 (
646 mesh_.faceCentres()[facei],
647 masterPointi,
648 -1,
649 true
650 );
651 }
652 }
653 }
654
655
656 // Per face, per point (faceCentre) or triangle (faceDiag) the (existing
657 // or added) cell on either side
658 faceOwnerCells_.setSize(mesh_.nFaces());
659 faceNeighbourCells_.setSize(mesh_.nFaces());
660
661 if (decomposeType == FACE_CENTRE_TRIS)
662 {
663 forAll(faceOwnerCells_, facei)
664 {
665 if (decomposeFace[facei])
666 {
667 const face& f = mesh_.faces()[facei];
668 faceOwnerCells_[facei].setSize(f.size(), -1);
669 faceNeighbourCells_[facei].setSize(f.size(), -1);
670 }
671 else
672 {
673 faceOwnerCells_[facei].setSize(1, -1);
674 faceNeighbourCells_[facei].setSize(1, -1);
675 }
676 }
677 }
678 else if (decomposeType == FACE_DIAG_TRIS)
679 {
680 // Force construction of diagonal decomposition
681 (void)mesh_.tetBasePtIs();
682
683 forAll(faceOwnerCells_, facei)
684 {
685 if (decomposeFace[facei])
686 {
687 const face& f = mesh_.faces()[facei];
688 faceOwnerCells_[facei].setSize(f.size()-2, -1);
689 faceNeighbourCells_[facei].setSize(f.size()-2, -1);
690 }
691 else
692 {
693 faceOwnerCells_[facei].setSize(1, -1);
694 faceNeighbourCells_[facei].setSize(1, -1);
695 }
696 }
697 }
698 else if (decomposeType == FACE_DIAG_QUADS)
699 {
700 // Note: sizing according to the number of sub-faces, not according to
701 // f.size(). Plus: less storage. Min: much harder to look up
702 // cell corresponding to face+edge
703
704 // Do coupled faces first - do not use local face triangulation
705 // but use coupled version instead
706 const auto& pbm = mesh_.boundaryMesh();
707 for (const auto& pp : pbm)
708 {
709 if (pp.coupled())
710 {
711 forAll(pp, i)
712 {
713 const label facei = pp.start()+i;
714 if (decomposeFace[facei])
715 {
716 const label bFacei = pp.offset()+i;
717 const label n =
718 boundaryQuads[bFacei].size()
719 + boundaryTris[bFacei].size();
720
721 faceOwnerCells_[facei].setSize(n, -1);
722 faceNeighbourCells_[facei].setSize(n, -1);
723 }
724 }
725 }
726 }
727
728 // Do all other faces
729 forAll(faceOwnerCells_, facei)
730 {
731 if (faceOwnerCells_[facei].empty())
732 {
733 if (decomposeFace[facei])
734 {
735 const face& f = mesh_.faces()[facei];
736
737 // Convention: quads first, followed by triangles
738
739 label nTris = 0;
740 label nQuads = 0;
741 const label nSubFaces = f.nTrianglesQuads
742 (
743 mesh_.points(),
744 nTris,
745 nQuads
746 );
747
748 faceOwnerCells_[facei].setSize(nSubFaces, -1);
749 faceNeighbourCells_[facei].setSize(nSubFaces, -1);
750 }
751 else
752 {
753 faceOwnerCells_[facei].setSize(1, -1);
754 faceNeighbourCells_[facei].setSize(1, -1);
755 }
756 }
757 }
758 }
759 else
760 {
761 forAll(faceOwnerCells_, facei)
762 {
763 faceOwnerCells_[facei].setSize(1, -1);
764 faceNeighbourCells_[facei].setSize(1, -1);
765 }
766 }
767
768
769 // Add internal cells. Note: done in same order as pyramid triangle
770 // creation later to maintain same ordering.
771 forAll(mesh_.cells(), celli)
772 {
773 const cell& cFaces = mesh_.cells()[celli];
774
775 // Whether cell has already been modified (all other cells get added)
776 bool modifiedCell = false;
777
778 forAll(cFaces, cFacei)
779 {
780 label facei = cFaces[cFacei];
781
782 // Get reference to either owner or neighbour
783 labelList& added =
784 (
785 (mesh_.faceOwner()[facei] == celli)
786 ? faceOwnerCells_[facei]
787 : faceNeighbourCells_[facei]
788 );
789
790 if (decomposeCell[celli])
791 {
792 if (decomposeType == FACE_CENTRE_TRIS)
793 {
794 forAll(added, i)
795 {
796 if (!modifiedCell)
797 {
798 // Reuse cell itself
799 added[i] = celli;
800 modifiedCell = true;
801 }
802 else
803 {
804 added[i] = meshMod.addCell
805 (
806 -1, // masterPoint
807 -1, // masterEdge
808 -1, // masterFace
809 celli, // masterCell
810 mesh_.cellZones().whichZone(celli)
811 );
812 }
813 }
814 }
815 else if
816 (
817 decomposeType == FACE_DIAG_TRIS
818 || decomposeType == FACE_DIAG_QUADS
819 )
820 {
821 forAll(added, subi)
822 {
823 if (!modifiedCell)
824 {
825 // Reuse cell itself
826 added[subi] = celli;
827 modifiedCell = true;
828 }
829 else
830 {
831 added[subi] = meshMod.addCell
832 (
833 -1, // masterPoint
834 -1, // masterEdge
835 -1, // masterFace
836 celli, // masterCell
837 mesh_.cellZones().whichZone(celli)
838 );
839 }
840 }
841 }
842 else // if (decomposeType == PYRAMID)
843 {
844 // Pyramidal decomposition.
845 // Assign same cell to all slots
846
847 if (added.size() != 1)
848 {
850 << "For cell:" << celli
851 << " at:" << mesh_.cellCentres()[celli]
852 << " face:" << facei
853 << " at:" << mesh_.faceCentres()[facei]
854 << " added cells:" << added
855 << exit(FatalError);
856 }
857
858
859 if (!modifiedCell)
860 {
861 // Reuse cell itself
862 added = celli;
863 modifiedCell = true;
864 }
865 else
866 {
867 added = meshMod.addCell
868 (
869 -1, // masterPoint
870 -1, // masterEdge
871 -1, // masterFace
872 celli, // masterCell
873 mesh_.cellZones().whichZone(celli)
874 );
875 }
876 }
877 }
878 else
879 {
880 // All vertices/triangles address to original cell
881 added = celli;
882 }
883 }
884 }
885
886
887
888 // Add split faces
889 face triangle(3);
890
891 forAll(mesh_.faces(), facei)
892 {
893 label own = mesh_.faceOwner()[facei];
894 const auto& addedOwn = faceOwnerCells_[facei];
895 const auto& addedNei = faceNeighbourCells_[facei];
896 const face& f = mesh_.faces()[facei];
897 //const point& fc = mesh_.faceCentres()[facei];
898
899 label nei = -1;
900 label patchi = -1;
901 if (facei >= mesh_.nInternalFaces())
902 {
903 patchi = mesh_.boundaryMesh().whichPatch(facei);
904 }
905 else
906 {
907 nei = mesh_.faceNeighbour()[facei];
908 }
909
910 label zonei = mesh_.faceZones().whichZone(facei);
911 bool zoneFlip = false;
912 if (zonei != -1)
913 {
914 const faceZone& fz = mesh_.faceZones()[zonei];
915 zoneFlip = fz.flipMap()[fz.whichFace(facei)];
916 }
917
918 if (decomposeType == FACE_CENTRE_TRIS && decomposeFace[facei])
919 {
920 forAll(f, fp)
921 {
922 // 1. Front triangle (decomposition of face itself)
923 // (between owner and neighbour cell)
924 {
925 triangle[0] = f[fp];
926 triangle[1] = f[f.fcIndex(fp)];
927 triangle[2] = faceToPoint_[facei];
928
929 const label newOwn
930 (
931 addedOwn.size() == 1
932 ? addedOwn[0]
933 : addedOwn[fp]
934 );
935 const label newNei
936 (
937 addedNei.size() == 1
938 ? addedNei[0]
939 : addedNei[fp]
940 );
941
942 if (fp == 0)
943 {
944 modifyFace
945 (
946 meshMod,
947 triangle,
948 facei,
949 newOwn,
950 newNei,
951 patchi,
952 zonei,
953 zoneFlip
954 );
955 }
956 else
957 {
958 addFace
959 (
960 meshMod,
961 triangle,
962 facei,
963 newOwn,
964 newNei,
965 -1, //point
966 -1, //edge
967 facei, //face
968 patchi,
969 zonei,
970 zoneFlip
971 );
972 }
973 }
974
975
976 // 2. Within owner cell - to cell centre
977 if (decomposeCell[own])
978 {
979 label newOwn = addedOwn[f.rcIndex(fp)];
980 label newNei = addedOwn[fp];
981
982 triangle[0] = f[fp];
983 triangle[1] = cellToPoint_[own];
984 triangle[2] = faceToPoint_[facei];
985
986 addFace
987 (
988 meshMod,
989 triangle,
990 facei,
991 newOwn,
992 newNei,
993 f[fp], //point
994 -1, //edge
995 -1, //face
996 -1, //patchi
997 -1, //zone
998 false
999 );
1000 }
1001 // 2b. Within neighbour cell - to cell centre
1002 if (facei < mesh_.nInternalFaces() && decomposeCell[nei])
1003 {
1004 label newOwn = addedNei[f.rcIndex(fp)];
1005 label newNei = addedNei[fp];
1006
1007 triangle[0] = f[fp];
1008 triangle[1] = faceToPoint_[facei];
1009 triangle[2] = cellToPoint_[nei];
1010
1011 addFace
1012 (
1013 meshMod,
1014 triangle,
1015 facei,
1016 newOwn,
1017 newNei,
1018 f[fp], //point
1019 -1, //edge
1020 -1, //face
1021 -1, //patchi
1022 -1, //zone
1023 false
1024 );
1025 }
1026 }
1027 }
1028 else if (decomposeType == FACE_DIAG_TRIS && decomposeFace[facei])
1029 {
1030 label fp0 = max(mesh_.tetBasePtIs()[facei], 0);
1031 label fp = f.fcIndex(fp0);
1032
1033 for (label trii = 0; trii < f.size()-2; trii++)
1034 {
1035 label nextTri = trii+1;
1036 if (nextTri >= f.size()-2)
1037 {
1038 nextTri -= f.size()-2;
1039 }
1040 label nextFp = f.fcIndex(fp);
1041
1042
1043 // Triangle trii consisiting of f[fp0], f[fp], f[nextFp]
1044
1045
1046 // 1. Front triangle (decomposition of face itself)
1047 // (between owner and neighbour cell)
1048 {
1049 const label newOwn
1050 (
1051 addedOwn.size() == 1
1052 ? addedOwn[0]
1053 : addedOwn[trii]
1054 );
1055 const label newNei
1056 (
1057 addedNei.size() == 1
1058 ? addedNei[0]
1059 : addedNei[trii]
1060 );
1061
1062 triangle[0] = f[fp0];
1063 triangle[1] = f[fp];
1064 triangle[2] = f[nextFp];
1065
1066 if (trii == 0)
1067 {
1068 modifyFace
1069 (
1070 meshMod,
1071 triangle,
1072 facei,
1073 newOwn,
1074 newNei,
1075 patchi,
1076 zonei,
1077 zoneFlip
1078 );
1079 }
1080 else
1081 {
1082 addFace
1083 (
1084 meshMod,
1085 triangle,
1086 facei,
1087 newOwn,
1088 newNei,
1089 -1, //point
1090 -1, //edge
1091 facei, //face
1092 patchi,
1093 zonei,
1094 zoneFlip
1095 );
1096 }
1097 }
1098
1099
1100 // 2. Within owner cell - diagonal to cell centre
1101 if (trii < f.size()-3)
1102 {
1103 if (decomposeCell[own])
1104 {
1105 label newOwn = addedOwn[trii];
1106 label newNei = addedOwn[nextTri];
1107
1108 triangle[0] = f[fp0];
1109 triangle[1] = f[nextFp];
1110 triangle[2] = cellToPoint_[own];
1111
1112 addFace
1113 (
1114 meshMod,
1115 triangle,
1116 facei,
1117 newOwn,
1118 newNei,
1119 f[fp], //point
1120 -1, //edge
1121 -1, //face
1122 -1, //patchi
1123 -1, //zone
1124 false
1125 );
1126 }
1127
1128 // 2b. Within neighbour cell - to cell centre
1129 if (facei < mesh_.nInternalFaces() && decomposeCell[nei])
1130 {
1131 label newOwn = addedNei[trii];
1132 label newNei = addedNei[nextTri];
1133
1134 triangle[0] = f[nextFp];
1135 triangle[1] = f[fp0];
1136 triangle[2] = cellToPoint_[nei];
1137
1138 addFace
1139 (
1140 meshMod,
1141 triangle,
1142 facei,
1143 newOwn,
1144 newNei,
1145 f[fp], //point
1146 -1, //edge
1147 -1, //face
1148 -1, //patchi
1149 -1, //zone
1150 false
1151 );
1152 }
1153 }
1154
1155
1156 fp = nextFp;
1157 }
1158 }
1159 else if (decomposeType == FACE_DIAG_QUADS && decomposeFace[facei])
1160 {
1161 // Decompose face into subFaces (quads followed by any triangles)
1162 splitFace
1163 (
1164 boundaryQuads,
1165 boundaryTris,
1166 facei,
1167 patchi,
1168 quadFaces, // work space
1169 triFaces, // work space
1170 subFaces
1171 );
1172
1173 EdgeMap<labelPair> edgeFaces(subFaces.size());
1174
1175 forAll(subFaces, subFacei)
1176 {
1177 const face& subF = subFaces[subFacei];
1178
1179 // 1. Front quad (decomposition of face itself)
1180 // (between owner and neighbour cell)
1181 {
1182 const label newOwn
1183 (
1184 addedOwn.size() == 1
1185 ? addedOwn[0]
1186 : addedOwn[subFacei]
1187 );
1188 const label newNei
1189 (
1190 addedNei.size() == 1
1191 ? addedNei[0]
1192 : addedNei[subFacei]
1193 );
1194
1195 if (subFacei == 0)
1196 {
1197 modifyFace
1198 (
1199 meshMod,
1200 subF,
1201 facei,
1202 newOwn,
1203 newNei,
1204 patchi,
1205 zonei,
1206 zoneFlip
1207 );
1208 }
1209 else
1210 {
1211 addFace
1212 (
1213 meshMod,
1214 subF,
1215 facei,
1216 newOwn,
1217 newNei,
1218 -1, //point
1219 -1, //edge
1220 facei, //face
1221 patchi,
1222 zonei,
1223 zoneFlip
1224 );
1225 }
1226 }
1227
1228 // Populate edge-faces (note: in order of increasing subFacei
1229 // and hence in order of added cells)
1230 forAll(subF, fp)
1231 {
1232 const edge e(subF.edge(fp));
1233 auto eFnd = edgeFaces.find(e);
1234 if (!eFnd)
1235 {
1236 edgeFaces.insert(e, labelPair(subFacei, -1));
1237 }
1238 else
1239 {
1240 auto& eFaces = eFnd();
1241 if (eFaces[1] != -1)
1242 {
1243 FatalErrorInFunction << "edgeFaces:" << edgeFaces
1244 << exit(FatalError);
1245 }
1246 eFaces[1] = subFacei;
1247 }
1248 }
1249 }
1250
1251 // Get diagonals
1252 forAllConstIters(edgeFaces, iter)
1253 {
1254 const auto& e = iter.key();
1255 const auto& eFaces = iter();
1256
1257 if (eFaces.find(-1) != -1)
1258 {
1259 // Open edge
1260 //Pout<< "for face:" << facei
1261 // << " ignoring open edge:" << e << endl;
1262 continue;
1263 }
1264
1265 if (decomposeCell[own])
1266 {
1267 const label newOwn(addedOwn[eFaces[0]]);
1268 const label newNei(addedOwn[eFaces[1]]);
1269
1270 if (newNei < newOwn) FatalErrorInFunction << "problem"
1271 << exit(FatalError);
1272
1273 // Point out of owner side
1274 triangle[0] = e[1];
1275 triangle[1] = e[0];
1276 triangle[2] = cellToPoint_[own];
1277
1278 addFace
1279 (
1280 meshMod,
1281 triangle,
1282 facei,
1283 newOwn,
1284 newNei,
1285 e[0], //point ? or edge ? or face ?
1286 -1, //edge
1287 -1, //face
1288 -1, //patchi
1289 -1, //zone
1290 false
1291 );
1292 }
1293
1294 // 2b. Within neighbour cell - to cell centre
1295 if
1296 (
1297 facei < mesh_.nInternalFaces()
1298 && decomposeCell[nei]
1299 )
1300 {
1301 const label newOwn(addedNei[eFaces[0]]);
1302 const label newNei(addedNei[eFaces[1]]);
1303
1304 if (newNei < newOwn) FatalErrorInFunction << "problem"
1305 << exit(FatalError);
1306
1307 triangle[0] = e[0];
1308 triangle[1] = e[1];
1309 triangle[2] = cellToPoint_[nei];
1310
1311 addFace
1312 (
1313 meshMod,
1314 triangle,
1315 facei,
1316 newOwn,
1317 newNei,
1318 e[0], //point ? or edge ? or face ?
1319 -1, //edge
1320 -1, //face
1321 -1, //patchi
1322 -1, //zone
1323 false
1324 );
1325 }
1326 }
1327 }
1328 else
1329 {
1330 // No decomposition. Use zero'th slot.
1331 modifyFace
1332 (
1333 meshMod,
1334 f,
1335 facei,
1336 addedOwn[0],
1337 addedNei[0],
1338 patchi,
1339 zonei,
1340 zoneFlip
1341 );
1342 }
1343 }
1344
1345
1346 // Add triangles to the cell centre for all edges to form the pyramids
1347 EdgeMap<label> edgeToFace;
1348
1349 forAll(mesh_.cells(), celli)
1350 {
1351 if (decomposeCell[celli])
1352 {
1353 const cell& cFaces = mesh_.cells()[celli];
1354
1355 edgeToFace.clear();
1356
1357 forAll(cFaces, cFacei)
1358 {
1359 const label facei = cFaces[cFacei];
1360 const label patchi = mesh_.boundaryMesh().whichPatch(facei);
1361 const face& f = mesh_.faces()[facei];
1362
1363 // Loop over edges of face. Avoid constructing faceEdges for
1364 // whole mesh.
1365 forAll(f, fp)
1366 {
1367 label p0 = f[fp];
1368 label p1 = f[f.fcIndex(fp)];
1369 const edge e(p0, p1);
1370
1371 EdgeMap<label>::const_iterator edgeFnd = edgeToFace.find(e);
1372 if (edgeFnd == edgeToFace.end())
1373 {
1374 edgeToFace.insert(e, facei);
1375 }
1376 else
1377 {
1378 // Found the other face on the edge.
1379 label otherFacei = edgeFnd();
1380 const face& otherF = mesh_.faces()[otherFacei];
1381
1382 // Found the other face on the edge. Note that since
1383 // we are looping in the same order the tets added for
1384 // otherFacei will be before those of facei
1385
1386 label otherFp = otherF.find(p0);
1387 if (otherF.nextLabel(otherFp) == p1)
1388 {
1389 // ok. otherFp is first vertex of edge.
1390 }
1391 else if (otherF.prevLabel(otherFp) == p1)
1392 {
1393 otherFp = otherF.rcIndex(otherFp);
1394 }
1395 else
1396 {
1398 << "problem." << abort(FatalError);
1399 }
1400
1401
1402 // Triangle from edge to cell centre
1403 if (mesh_.faceOwner()[facei] == celli)
1404 {
1405 triangle[0] = p0;
1406 triangle[1] = p1;
1407 triangle[2] = cellToPoint_[celli];
1408 }
1409 else
1410 {
1411 triangle[0] = p1;
1412 triangle[1] = p0;
1413 triangle[2] = cellToPoint_[celli];
1414 }
1415
1416 // Determine tets on either side
1417
1418 const auto& thisCells =
1419 (
1420 (mesh_.faceOwner()[facei] == celli)
1421 ? faceOwnerCells_[facei]
1422 : faceNeighbourCells_[facei]
1423 );
1424 const auto& otherCells =
1425 (
1426 (mesh_.faceOwner()[otherFacei] == celli)
1427 ? faceOwnerCells_[otherFacei]
1428 : faceNeighbourCells_[otherFacei]
1429 );
1430
1431 label thisTet = -1;
1432 label otherTet = -1;
1433
1434 if (decomposeType == FACE_CENTRE_TRIS)
1435 {
1436 if (thisCells.size() == 1)
1437 {
1438 thisTet = thisCells[0];
1439 }
1440 else
1441 {
1442 thisTet = thisCells[fp];
1443 }
1444
1445 if (otherCells.size() == 1)
1446 {
1447 otherTet = otherCells[0];
1448 }
1449 else
1450 {
1451 otherTet = otherCells[otherFp];
1452 }
1453 }
1454 else if (decomposeType == FACE_DIAG_TRIS)
1455 {
1456 if (thisCells.size() == 1)
1457 {
1458 thisTet = thisCells[0];
1459 }
1460 else
1461 {
1462 thisTet = thisCells[triIndex(facei, fp)];
1463 }
1464
1465 if (otherCells.size() == 1)
1466 {
1467 otherTet = otherCells[0];
1468 }
1469 else
1470 {
1471 otherTet =
1472 otherCells[triIndex(otherFacei, otherFp)];
1473 }
1474 }
1475 else if (decomposeType == FACE_DIAG_QUADS)
1476 {
1477 // Get added cell on facei
1478 if (thisCells.size() == 1)
1479 {
1480 thisTet = thisCells[0];
1481 }
1482 else
1483 {
1484 // We have not stored the face-split so
1485 // don't know yet the correspondence. Maybe
1486 // store to avoid re-doing split?
1487 splitFace
1488 (
1489 boundaryQuads,
1490 boundaryTris,
1491 facei,
1492 patchi,
1493 quadFaces, // work space
1494 triFaces, // work space
1495 subFaces
1496 );
1497
1498 // Find the edge. Should be only on one of the
1499 // subfaces.
1500 forAll(subFaces, subFacei)
1501 {
1502 const auto& subF = subFaces[subFacei];
1503 if (subF.find(e) != -1)
1504 {
1505 thisTet = thisCells[subFacei];
1506 break;
1507 }
1508 }
1509 if (thisTet == -1)
1510 {
1512 << "For cell:" << celli
1513 << " at:" << mesh_.cellCentres()[celli]
1514 << " patch:" << patchi
1515 << " face:" << facei
1516 << " at:" << mesh_.faceCentres()[facei]
1517 << " did not find edge:" << e
1518 << " at:" << e.line(mesh_.points())
1519 << " in OWNER-side decomposed faces:"
1520 << flatOutput(subFaces)
1521 << exit(FatalError);
1522 }
1523 }
1524
1525
1526 // Get added cell on otherFacei
1527 if (otherCells.size() == 1)
1528 {
1529 otherTet = otherCells[0];
1530 }
1531 else
1532 {
1533 splitFace
1534 (
1535 boundaryQuads,
1536 boundaryTris,
1537 otherFacei,
1538 mesh_.boundaryMesh().whichPatch(otherFacei),
1539 quadFaces, // work space
1540 triFaces, // work space
1541 subFaces
1542 );
1543
1544 // Find the edge. Should be only on one of the
1545 // subfaces.
1546 forAll(subFaces, subFacei)
1547 {
1548 const auto& subF = subFaces[subFacei];
1549 if (subF.find(e) != -1)
1550 {
1551 otherTet = otherCells[subFacei];
1552 break;
1553 }
1554 }
1555 if (otherTet == -1)
1556 {
1558 << "For cell:" << celli
1559 << " at:" << mesh_.cellCentres()[celli]
1560 << " patch:"
1561 << mesh_.boundaryMesh().whichPatch(otherFacei)
1562 << " face:" << otherFacei
1563 << " at:"
1564 << mesh_.faceCentres()[otherFacei]
1565 << " did not find edge:" << e
1566 << " at:" << e.line(mesh_.points())
1567 << " in decomposed faces:"
1568 << flatOutput(subFaces)
1569 << exit(FatalError);
1570 }
1571 }
1572 }
1573 else
1574 {
1575 thisTet = thisCells[0];
1576 otherTet = otherCells[0];
1577 }
1578
1579 addFace
1580 (
1581 meshMod,
1582 triangle,
1583 facei,
1584 otherTet,
1585 thisTet,
1586 -1, //masterPoint
1587 -1, //fEdges[fp], //masterEdge
1588 facei, //masterFace
1589 -1, //patchi
1590 -1, //zone
1591 false
1592 );
1593 }
1595 }
1596 }
1597 }
1598}
1599
1600
1602{
1603 inplaceRenumber(map.reversePointMap(), cellToPoint_);
1604 inplaceRenumber(map.reversePointMap(), faceToPoint_);
1605
1606 forAll(faceOwnerCells_, facei)
1607 {
1608 inplaceRenumber(map.reverseCellMap(), faceOwnerCells_[facei]);
1609 }
1610 forAll(faceNeighbourCells_, facei)
1611 {
1612 inplaceRenumber(map.reverseCellMap(), faceNeighbourCells_[facei]);
1613 }
1614}
1615
1616
1617// ************************************************************************* //
label n
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
Map from edge (expressed as its endpoints) to value. Hashing (and ==) on an edge is symmetric.
Definition edgeHashes.H:59
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition HashTableI.H:86
void clear()
Remove all entries from table.
Definition HashTable.C:742
iterator end() noexcept
iterator to signal the end (for any HashTable)
label size() const noexcept
The number of elements in the list.
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
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
void setSize(label n)
Alias for resize().
Definition List.H:536
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
label rcIndex(const label i) const noexcept
The reverse circular index. The previous index in the list which returns to the last at the beginning...
Definition UListI.H:106
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition UListI.H:99
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
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 subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
label whichFace(const label meshFaceID) const
The local index of the given mesh face, -1 if not in the zone.
Definition faceZone.H:394
const boolList & flipMap() const noexcept
Return face flip map.
Definition faceZone.H:389
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Foam::edge edge(const label edgei) const
Return i-th face edge (forward walk order).
Definition faceI.H:158
label find(const Foam::edge &e) const
Find the edge within the face.
Definition face.C:791
label nextLabel(const label i) const
Next vertex on face.
Definition faceI.H:208
label prevLabel(const label i) const
Previous vertex on face.
Definition faceI.H:214
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reverseCellMap() const noexcept
Reverse cell map.
const labelList & reversePointMap() const noexcept
Reverse point map.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Direct mesh changes based on v1.3 polyTopoChange syntax.
label addPoint(const point &pt, const label masterPointID, const label zoneID, const bool inCell)
Add point. Return new point label.
label addCell(const label masterPointID, const label masterEdgeID, const label masterFaceID, const label masterCellID, const label zoneID)
Add cell. Return new cell label.
static void syncBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const TransformOp &top, const bool parRun=UPstream::parRun())
Synchronize values on boundary faces only.
static void swapFaceList(const polyMesh &mesh, UList< T > &faceValues, const bool parRun=UPstream::parRun())
Swap coupled face values. Uses eqOp.
Definition syncTools.H:567
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const bool parRun=UPstream::parRun())
Swap coupled boundary face values. Uses eqOp.
Definition syncTools.H:524
Decomposes polyMesh into tets (or pyramids).
static const Enum< decompositionType > decompositionTypeNames
void setRefinement(const decompositionType decomposeType, const bitSet &decomposeCell, polyTopoChange &meshMod)
Insert all changes into meshMod to convert the polyMesh into.
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
A triangle primitive used to calculate face normals and swept volumes. Uses referred points.
Definition triangle.H:234
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const volScalarField & p0
Definition EEqn.H:36
dynamicFvMesh & mesh
Dummy transform to be used with syncTools.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
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< label > labelList
A List of labels.
Definition List.H:62
List< face > faceList
List of faces.
Definition faceListFwd.H:41
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
errorManip< error > abort(error &err)
Definition errorManip.H:139
List< bool > boolList
A List of bools.
Definition List.H:60
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
labelList f(nPoints)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition stdFoam.H:315
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235