meshCutAndRemove.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) 2019 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 "meshCutAndRemove.H"
30#include "polyMesh.H"
31#include "polyTopoChange.H"
32#include "polyAddFace.H"
33#include "polyAddPoint.H"
34#include "polyRemovePoint.H"
35#include "polyRemoveFace.H"
36#include "polyModifyFace.H"
37#include "cellCuts.H"
38#include "mapPolyMesh.H"
39#include "meshTools.H"
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43namespace Foam
44{
46}
47
48
49// * * * * * * * * * * * * * Private Static Functions * * * * * * * * * * * //
50
51// Returns -1 or index in elems1 of first shared element.
52Foam::label Foam::meshCutAndRemove::firstCommon
53(
54 const labelList& elems1,
55 const labelList& elems2
56)
57{
58 forAll(elems1, elemI)
59 {
60 label index1 = elems2.find(elems1[elemI]);
61
62 if (index1 != -1)
63 {
64 return index1;
65 }
66 }
67 return -1;
68}
69
70
71// Check if twoCuts at two consecutive position in cuts.
72bool Foam::meshCutAndRemove::isIn
73(
74 const edge& twoCuts,
75 const labelList& cuts
76)
77{
78 label index = cuts.find(twoCuts[0]);
79
80 if (index == -1)
81 {
82 return false;
83 }
84
85 return
86 (
87 cuts[cuts.fcIndex(index)] == twoCuts[1]
88 || cuts[cuts.rcIndex(index)] == twoCuts[1]
89 );
90}
91
92
93// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
94
95Foam::label Foam::meshCutAndRemove::findCutCell
96(
97 const cellCuts& cuts,
98 const labelList& cellLabels
99) const
100{
101 forAll(cellLabels, labelI)
102 {
103 label celli = cellLabels[labelI];
104
105 if (cuts.cellLoops()[celli].size())
106 {
107 return celli;
108 }
109 }
110 return -1;
111}
112
113
114Foam::label Foam::meshCutAndRemove::findInternalFacePoint
115(
117) const
118{
119 forAll(pointLabels, labelI)
120 {
121 label pointi = pointLabels[labelI];
122
123 const labelList& pFaces = mesh().pointFaces()[pointi];
124
125 forAll(pFaces, pFacei)
126 {
127 label facei = pFaces[pFacei];
128
129 if (mesh().isInternalFace(facei))
130 {
131 return pointi;
132 }
133 }
134 }
135
136 if (pointLabels.empty())
137 {
139 << "Empty pointLabels" << abort(FatalError);
140 }
141
142 return -1;
143}
144
145
146Foam::label Foam::meshCutAndRemove::findPatchFacePoint
147(
148 const face& f,
149 const label exposedPatchi
150) const
151{
152 const labelListList& pointFaces = mesh().pointFaces();
153 const polyBoundaryMesh& patches = mesh().boundaryMesh();
154
155 forAll(f, fp)
156 {
157 label pointi = f[fp];
158
159 if (pointi < mesh().nPoints())
160 {
161 const labelList& pFaces = pointFaces[pointi];
162
163 forAll(pFaces, i)
164 {
165 if (patches.whichPatch(pFaces[i]) == exposedPatchi)
166 {
167 return pointi;
168 }
169 }
170 }
171 }
172 return -1;
173}
174
175
177(
178 const cellCuts& cuts,
179 const label exposedPatchi,
180 const label facei,
181 label& own,
182 label& nei,
183 label& patchID
184) const
185{
186 const labelListList& anchorPts = cuts.cellAnchorPoints();
187 const labelListList& cellLoops = cuts.cellLoops();
188
189 const face& f = mesh().faces()[facei];
190
191 own = mesh().faceOwner()[facei];
192
193 if (cellLoops[own].size() && firstCommon(f, anchorPts[own]) == -1)
194 {
195 // owner has been split and this is the removed part.
196 own = -1;
197 }
198
199 nei = -1;
200
201 if (mesh().isInternalFace(facei))
202 {
203 nei = mesh().faceNeighbour()[facei];
204
205 if (cellLoops[nei].size() && firstCommon(f, anchorPts[nei]) == -1)
206 {
207 nei = -1;
208 }
209 }
210
211 patchID = mesh().boundaryMesh().whichPatch(facei);
212
213 if (patchID == -1 && (own == -1 || nei == -1))
214 {
215 // Face was internal but becomes external
216 patchID = exposedPatchi;
217 }
218}
219
220
221void Foam::meshCutAndRemove::getZoneInfo
222(
223 const label facei,
224 label& zoneID,
225 bool& zoneFlip
226) const
227{
228 zoneID = mesh().faceZones().whichZone(facei);
229
230 zoneFlip = false;
231
232 if (zoneID >= 0)
233 {
234 const faceZone& fZone = mesh().faceZones()[zoneID];
235
236 zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
237 }
238}
239
240
241void Foam::meshCutAndRemove::addFace
242(
243 polyTopoChange& meshMod,
244 const label facei,
245 const label masterPointi,
246 const face& newFace,
247 const label own,
248 const label nei,
249 const label patchID
250)
251{
252 label zoneID;
253 bool zoneFlip;
254
255 getZoneInfo(facei, zoneID, zoneFlip);
256
257 if ((nei == -1) || (own != -1 && own < nei))
258 {
259 // Ordering ok.
260 if (debug & 2)
261 {
262 Pout<< "Adding face " << newFace
263 << " with new owner:" << own
264 << " with new neighbour:" << nei
265 << " patchID:" << patchID
266 << " anchor:" << masterPointi
267 << " zoneID:" << zoneID
268 << " zoneFlip:" << zoneFlip
269 << endl;
270 }
271
272 meshMod.setAction
273 (
274 polyAddFace
275 (
276 newFace, // face
277 own, // owner
278 nei, // neighbour
279 masterPointi, // master point
280 -1, // master edge
281 -1, // master face for addition
282 false, // flux flip
283 patchID, // patch for face
284 zoneID, // zone for face
285 zoneFlip // face zone flip
286 )
287 );
288 }
289 else
290 {
291 // Reverse owner/neighbour
292 if (debug & 2)
293 {
294 Pout<< "Adding (reversed) face " << newFace.reverseFace()
295 << " with new owner:" << nei
296 << " with new neighbour:" << own
297 << " patchID:" << patchID
298 << " anchor:" << masterPointi
299 << " zoneID:" << zoneID
300 << " zoneFlip:" << zoneFlip
301 << endl;
302 }
303
304 meshMod.setAction
305 (
306 polyAddFace
307 (
308 newFace.reverseFace(), // face
309 nei, // owner
310 own, // neighbour
311 masterPointi, // master point
312 -1, // master edge
313 -1, // master face for addition
314 false, // flux flip
315 patchID, // patch for face
316 zoneID, // zone for face
317 zoneFlip // face zone flip
318 )
319 );
320 }
321}
322
323
324// Modifies existing facei for either new owner/neighbour or new face points.
325void Foam::meshCutAndRemove::modFace
326(
327 polyTopoChange& meshMod,
328 const label facei,
329 const face& newFace,
330 const label own,
331 const label nei,
332 const label patchID
333)
334{
335 label zoneID;
336 bool zoneFlip;
337
338 getZoneInfo(facei, zoneID, zoneFlip);
339
340 if
341 (
342 (own != mesh().faceOwner()[facei])
343 || (
344 mesh().isInternalFace(facei)
345 && (nei != mesh().faceNeighbour()[facei])
346 )
347 || (newFace != mesh().faces()[facei])
348 )
349 {
350 if (debug & 2)
351 {
352 Pout<< "Modifying face " << facei
353 << " old vertices:" << mesh().faces()[facei]
354 << " new vertices:" << newFace
355 << " new owner:" << own
356 << " new neighbour:" << nei
357 << " new patch:" << patchID
358 << " new zoneID:" << zoneID
359 << " new zoneFlip:" << zoneFlip
360 << endl;
361 }
362
363 if ((nei == -1) || (own != -1 && own < nei))
364 {
365 meshMod.setAction
366 (
367 polyModifyFace
368 (
369 newFace, // modified face
370 facei, // label of face being modified
371 own, // owner
372 nei, // neighbour
373 false, // face flip
374 patchID, // patch for face
375 false, // remove from zone
376 zoneID, // zone for face
377 zoneFlip // face flip in zone
378 )
379 );
380 }
381 else
382 {
383 meshMod.setAction
384 (
385 polyModifyFace
386 (
387 newFace.reverseFace(), // modified face
388 facei, // label of face being modified
389 nei, // owner
390 own, // neighbour
391 false, // face flip
392 patchID, // patch for face
393 false, // remove from zone
394 zoneID, // zone for face
395 zoneFlip // face flip in zone
396 )
397 );
398 }
399 }
400}
401
402
403void Foam::meshCutAndRemove::copyFace
404(
405 const face& f,
406 const label startFp,
407 const label endFp,
408 face& newFace
409) const
410{
411 label fp = startFp;
412
413 label newFp = 0;
414
415 while (fp != endFp)
416 {
417 newFace[newFp++] = f[fp];
418
419 fp = (fp + 1) % f.size();
420 }
421 newFace[newFp] = f[fp];
422}
423
424
425// Actually split face in two along splitEdge v0, v1 (the two vertices in new
426// vertex numbering). Generates faces in same ordering
427// as original face. Replaces cutEdges by the points introduced on them
428// (addedPoints_).
429void Foam::meshCutAndRemove::splitFace
430(
431 const face& f,
432 const label v0,
433 const label v1,
434
435 face& f0,
436 face& f1
437) const
438{
439 // Check if we find any new vertex which is part of the splitEdge.
440 label startFp = f.find(v0);
441
442 if (startFp == -1)
443 {
445 << "Cannot find vertex (new numbering) " << v0
446 << " on face " << f
447 << abort(FatalError);
448 }
449
450 label endFp = f.find(v1);
451
452 if (endFp == -1)
453 {
455 << "Cannot find vertex (new numbering) " << v1
456 << " on face " << f
457 << abort(FatalError);
458 }
459
460
461 f0.setSize((endFp + 1 + f.size() - startFp) % f.size());
462 f1.setSize(f.size() - f0.size() + 2);
463
464 copyFace(f, startFp, endFp, f0);
465 copyFace(f, endFp, startFp, f1);
466}
467
468
469Foam::face Foam::meshCutAndRemove::addEdgeCutsToFace(const label facei) const
470{
471 const face& f = mesh().faces()[facei];
472
473 face newFace(2 * f.size());
474
475 label newFp = 0;
476
477 forAll(f, fp)
478 {
479 // Duplicate face vertex.
480 newFace[newFp++] = f[fp];
481
482 // Check if edge has been cut.
483 label fp1 = f.fcIndex(fp);
484
485 EdgeMap<label>::const_iterator fnd =
486 addedPoints_.find(edge(f[fp], f[fp1]));
487
488 if (fnd.found())
489 {
490 // edge has been cut. Introduce new vertex.
491 newFace[newFp++] = fnd.val();
492 }
493 }
494
495 newFace.setSize(newFp);
496
497 return newFace;
498}
499
500
501// Walk loop (loop of cuts) across circumference of celli. Returns face in
502// new vertices.
503// Note: tricky bit is that it can use existing edges which have been split.
504Foam::face Foam::meshCutAndRemove::loopToFace
505(
506 const label celli,
507 const labelList& loop
508) const
509{
510 face newFace(2*loop.size());
511
512 label newFacei = 0;
513
514 forAll(loop, fp)
515 {
516 label cut = loop[fp];
517
518 if (isEdge(cut))
519 {
520 label edgeI = getEdge(cut);
521
522 const edge& e = mesh().edges()[edgeI];
523
524 label vertI = addedPoints_[e];
525
526 newFace[newFacei++] = vertI;
527 }
528 else
529 {
530 // cut is vertex.
531 label vertI = getVertex(cut);
532
533 newFace[newFacei++] = vertI;
534
535 label nextCut = loop[loop.fcIndex(fp)];
536
537 if (!isEdge(nextCut))
538 {
539 // From vertex to vertex -> cross cut only if no existing edge.
540 label nextVertI = getVertex(nextCut);
541
542 label edgeI = meshTools::findEdge(mesh(), vertI, nextVertI);
543
544 if (edgeI != -1)
545 {
546 // Existing edge. Insert split-edge point if any.
547 EdgeMap<label>::const_iterator fnd =
548 addedPoints_.find(mesh().edges()[edgeI]);
549
550 if (fnd.found())
551 {
552 newFace[newFacei++] = fnd.val();
553 }
554 }
555 }
556 }
557 }
558 newFace.setSize(newFacei);
559
560 return newFace;
561}
562
563
564// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
565
567:
569 addedFaces_(),
570 addedPoints_()
571{}
572
573
574// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
575
577(
578 const label exposedPatchi,
579 const cellCuts& cuts,
580 const labelList& cutPatch,
581 polyTopoChange& meshMod
582)
583{
584 // Clear and size maps here since mesh size will change.
585 addedFaces_.clear();
586 addedFaces_.resize(cuts.nLoops());
587
588 addedPoints_.clear();
589 addedPoints_.resize(cuts.nLoops());
590
591 if (cuts.nLoops() == 0)
592 {
593 return;
594 }
595
596 const labelListList& anchorPts = cuts.cellAnchorPoints();
597 const labelListList& cellLoops = cuts.cellLoops();
599
600 if (exposedPatchi < 0 || exposedPatchi >= patches.size())
601 {
603 << "Illegal exposed patch " << exposedPatchi
604 << abort(FatalError);
605 }
606
607
608 //
609 // Add new points along cut edges.
610 //
611
612 forAll(cuts.edgeIsCut(), edgeI)
613 {
614 if (cuts.edgeIsCut()[edgeI])
615 {
616 const edge& e = mesh().edges()[edgeI];
617
618 // Check if there is any cell using this edge.
619 if (debug && findCutCell(cuts, mesh().edgeCells()[edgeI]) == -1)
620 {
622 << "Problem: cut edge but none of the cells using it is\n"
623 << "edge:" << edgeI << " verts:" << e
624 << abort(FatalError);
625 }
626
627 // One of the edge end points should be master point of nbCelli.
628 label masterPointi = e.start();
629
630 const point& v0 = mesh().points()[e.start()];
631 const point& v1 = mesh().points()[e.end()];
632
633 scalar weight = cuts.edgeWeight()[edgeI];
634
635 point newPt = weight*v1 + (1.0-weight)*v0;
636
637 label addedPointi =
638 meshMod.setAction
639 (
641 (
642 newPt, // point
643 masterPointi, // master point
644 -1, // zone for point
645 true // supports a cell
646 )
647 );
648
649 // Store on (hash of) edge.
650 addedPoints_.insert(e, addedPointi);
651
652 if (debug & 2)
653 {
654 Pout<< "Added point " << addedPointi
655 << " to vertex "
656 << masterPointi << " of edge " << edgeI
657 << " vertices " << e << endl;
658 }
659 }
660 }
661
662
663 //
664 // Remove all points that will not be used anymore
665 //
666 {
667 boolList usedPoint(mesh().nPoints(), false);
668
669 forAll(cellLoops, celli)
670 {
671 const labelList& loop = cellLoops[celli];
672
673 if (loop.size())
674 {
675 // Cell is cut. Uses only anchor points and loop itself.
676 forAll(loop, fp)
677 {
678 label cut = loop[fp];
679
680 if (!isEdge(cut))
681 {
682 usedPoint[getVertex(cut)] = true;
683 }
684 }
685
686 const labelList& anchors = anchorPts[celli];
687
688 forAll(anchors, i)
689 {
690 usedPoint[anchors[i]] = true;
691 }
692 }
693 else
694 {
695 // Cell is not cut so use all its points
696 const labelList& cPoints = mesh().cellPoints()[celli];
697
698 forAll(cPoints, i)
699 {
700 usedPoint[cPoints[i]] = true;
701 }
702 }
703 }
704
705
706 // Check
707 const Map<edge>& faceSplitCut = cuts.faceSplitCut();
708
709 forAllConstIters(faceSplitCut, iter)
710 {
711 const edge& fCut = iter.val();
712
713 forAll(fCut, i)
714 {
715 label cut = fCut[i];
716
717 if (!isEdge(cut))
718 {
719 label pointi = getVertex(cut);
720
721 if (!usedPoint[pointi])
722 {
724 << "Problem: faceSplitCut not used by any loop"
725 << " or cell anchor point"
726 << "face:" << iter.key() << " point:" << pointi
727 << " coord:" << mesh().points()[pointi]
728 << abort(FatalError);
729 }
730 }
731 }
732 }
733
734 forAll(cuts.pointIsCut(), pointi)
735 {
736 if (cuts.pointIsCut()[pointi])
737 {
738 if (!usedPoint[pointi])
739 {
741 << "Problem: point is marked as cut but"
742 << " not used by any loop"
743 << " or cell anchor point"
744 << "point:" << pointi
745 << " coord:" << mesh().points()[pointi]
746 << abort(FatalError);
747 }
748 }
749 }
750
751
752 // Remove unused points.
753 forAll(usedPoint, pointi)
754 {
755 if (!usedPoint[pointi])
756 {
757 meshMod.setAction(polyRemovePoint(pointi));
758
759 if (debug & 2)
760 {
761 Pout<< "Removing unused point " << pointi << endl;
762 }
763 }
764 }
765 }
766
767
768 //
769 // For all cut cells add an internal or external face
770 //
771
772 forAll(cellLoops, celli)
773 {
774 const labelList& loop = cellLoops[celli];
775
776 if (loop.size())
777 {
778 if (cutPatch[celli] < 0 || cutPatch[celli] >= patches.size())
779 {
781 << "Illegal patch " << cutPatch[celli]
782 << " provided for cut cell " << celli
783 << abort(FatalError);
784 }
785
786 //
787 // Convert loop (=list of cuts) into proper face.
788 // cellCuts sets orientation is towards anchor side so reverse.
789 //
790 face newFace(loopToFace(celli, loop));
791
792 reverse(newFace);
793
794 // Pick any anchor point on cell
795 label masterPointi = findPatchFacePoint(newFace, exposedPatchi);
796
797 label addedFacei =
798 meshMod.setAction
799 (
801 (
802 newFace, // face
803 celli, // owner
804 -1, // neighbour
805 masterPointi, // master point
806 -1, // master edge
807 -1, // master face for addition
808 false, // flux flip
809 cutPatch[celli], // patch for face
810 -1, // zone for face
811 false // face zone flip
812 )
813 );
814
815 addedFaces_.insert(celli, addedFacei);
816
817 if (debug & 2)
818 {
819 Pout<< "Added splitting face " << newFace << " index:"
820 << addedFacei << " from masterPoint:" << masterPointi
821 << " to owner " << celli << " with anchors:"
822 << anchorPts[celli]
823 << " from Loop:";
824
825 // Gets edgeweights of loop
826 scalarField weights(loop.size());
827 forAll(loop, i)
828 {
829 label cut = loop[i];
830
831 weights[i] =
832 (
833 isEdge(cut)
834 ? cuts.edgeWeight()[getEdge(cut)]
835 : -GREAT
836 );
837 }
838 writeCuts(Pout, loop, weights);
839 Pout<< endl;
840 }
841 }
842 }
843
844
845 //
846 // Modify faces to use only anchorpoints and loop points
847 // (so throw away part without anchorpoints)
848 //
849
850
851 // Maintain whether face has been updated (for -split edges
852 // -new owner/neighbour)
853 boolList faceUptodate(mesh().nFaces(), false);
854
855 const Map<edge>& faceSplitCuts = cuts.faceSplitCut();
856
857 forAllConstIters(faceSplitCuts, iter)
858 {
859 const label facei = iter.key();
860
861 const edge& splitEdge = iter.val();
862
863 // Renumber face to include split edges.
864 face newFace(addEdgeCutsToFace(facei));
865
866 // Edge splitting the face. Convert edge to new vertex numbering.
867 label cut0 = splitEdge[0];
868
869 label v0;
870 if (isEdge(cut0))
871 {
872 label edgeI = getEdge(cut0);
873 v0 = addedPoints_[mesh().edges()[edgeI]];
874 }
875 else
876 {
877 v0 = getVertex(cut0);
878 }
879
880 label cut1 = splitEdge[1];
881 label v1;
882 if (isEdge(cut1))
883 {
884 label edgeI = getEdge(cut1);
885 v1 = addedPoints_[mesh().edges()[edgeI]];
886 }
887 else
888 {
889 v1 = getVertex(cut1);
890 }
891
892 // Split face along the elements of the splitEdge.
893 face f0, f1;
894 splitFace(newFace, v0, v1, f0, f1);
895
896 label own = mesh().faceOwner()[facei];
897
898 label nei = -1;
899
900 if (mesh().isInternalFace(facei))
901 {
902 nei = mesh().faceNeighbour()[facei];
903 }
904
905 if (debug & 2)
906 {
907 Pout<< "Split face " << mesh().faces()[facei]
908 << " own:" << own << " nei:" << nei
909 << " into f0:" << f0
910 << " and f1:" << f1 << endl;
911 }
912
913
914 // Check which cell using face uses anchorPoints (so is kept)
915 // and which one doesn't (gets removed)
916
917 // Bit tricky. We have to know whether this faceSplit splits owner/
918 // neighbour or both. Even if cell is cut we have to make sure this is
919 // the one that cuts it (this face cut might not be the one splitting
920 // the cell)
921 // The face f gets split into two parts, f0 and f1.
922 // Each of these can have a different owner and or neighbour.
923
924 const face& f = mesh().faces()[facei];
925
926 label f0Own = -1;
927 label f1Own = -1;
928
929 if (cellLoops[own].empty())
930 {
931 // Owner side is not split so keep both halves.
932 f0Own = own;
933 f1Own = own;
934 }
935 else if (isIn(splitEdge, cellLoops[own]))
936 {
937 // Owner is cut by this splitCut. See which of f0, f1 gets
938 // preserved and becomes owner, and which gets removed.
939 if (firstCommon(f0, anchorPts[own]) != -1)
940 {
941 // f0 preserved so f1 gets deleted
942 f0Own = own;
943 f1Own = -1;
944 }
945 else
946 {
947 f0Own = -1;
948 f1Own = own;
949 }
950 }
951 else
952 {
953 // Owner not cut by this splitCut but by another.
954 // Check on original face whether
955 // use anchorPts.
956 if (firstCommon(f, anchorPts[own]) != -1)
957 {
958 // both f0 and f1 owner side preserved
959 f0Own = own;
960 f1Own = own;
961 }
962 else
963 {
964 // both f0 and f1 owner side removed
965 f0Own = -1;
966 f1Own = -1;
967 }
968 }
969
970
971 label f0Nei = -1;
972 label f1Nei = -1;
973
974 if (nei != -1)
975 {
976 if (cellLoops[nei].empty())
977 {
978 f0Nei = nei;
979 f1Nei = nei;
980 }
981 else if (isIn(splitEdge, cellLoops[nei]))
982 {
983 // Neighbour is cut by this splitCut. So anchor part of it
984 // gets kept, non-anchor bit gets removed. See which of f0, f1
985 // connects to which part.
986
987 if (firstCommon(f0, anchorPts[nei]) != -1)
988 {
989 f0Nei = nei;
990 f1Nei = -1;
991 }
992 else
993 {
994 f0Nei = -1;
995 f1Nei = nei;
996 }
997 }
998 else
999 {
1000 // neighbour not cut by this splitCut. Check on original face
1001 // whether use anchorPts.
1002
1003 if (firstCommon(f, anchorPts[nei]) != -1)
1004 {
1005 f0Nei = nei;
1006 f1Nei = nei;
1007 }
1008 else
1009 {
1010 // both f0 and f1 on neighbour side removed
1011 f0Nei = -1;
1012 f1Nei = -1;
1013 }
1014 }
1015 }
1016
1017
1018 if (debug & 2)
1019 {
1020 Pout<< "f0 own:" << f0Own << " nei:" << f0Nei
1021 << " f1 own:" << f1Own << " nei:" << f1Nei
1022 << endl;
1023 }
1024
1025
1026 // If faces were internal but now become external set a patch.
1027 // If they were external already keep the patch.
1028 label patchID = patches.whichPatch(facei);
1029
1030 if (patchID == -1)
1031 {
1032 patchID = exposedPatchi;
1033 }
1034
1035
1036 // Do as much as possible by modifying facei. Delay any remove
1037 // face. Keep track of whether facei has been used.
1038
1039 bool modifiedFacei = false;
1040
1041 if (f0Own == -1)
1042 {
1043 if (f0Nei != -1)
1044 {
1045 // f0 becomes external face (note:modFace will reverse face)
1046 modFace(meshMod, facei, f0, f0Own, f0Nei, patchID);
1047 modifiedFacei = true;
1048 }
1049 }
1050 else
1051 {
1052 if (f0Nei == -1)
1053 {
1054 // f0 becomes external face
1055 modFace(meshMod, facei, f0, f0Own, f0Nei, patchID);
1056 modifiedFacei = true;
1057 }
1058 else
1059 {
1060 // f0 stays internal face.
1061 modFace(meshMod, facei, f0, f0Own, f0Nei, -1);
1062 modifiedFacei = true;
1063 }
1064 }
1065
1066
1067 // f1 is added face (if at all)
1068
1069 if (f1Own == -1)
1070 {
1071 if (f1Nei == -1)
1072 {
1073 // f1 not needed.
1074 }
1075 else
1076 {
1077 // f1 becomes external face (note:modFace will reverse face)
1078 if (!modifiedFacei)
1079 {
1080 modFace(meshMod, facei, f1, f1Own, f1Nei, patchID);
1081 modifiedFacei = true;
1082 }
1083 else
1084 {
1085 label masterPointi = findPatchFacePoint(f1, patchID);
1086
1087 addFace
1088 (
1089 meshMod,
1090 facei, // face for zone info
1091 masterPointi, // inflation point
1092 f1, // vertices of face
1093 f1Own,
1094 f1Nei,
1095 patchID // patch for new face
1096 );
1097 }
1098 }
1099 }
1100 else
1101 {
1102 if (f1Nei == -1)
1103 {
1104 // f1 becomes external face
1105 if (!modifiedFacei)
1106 {
1107 modFace(meshMod, facei, f1, f1Own, f1Nei, patchID);
1108 modifiedFacei = true;
1109 }
1110 else
1111 {
1112 label masterPointi = findPatchFacePoint(f1, patchID);
1113
1114 addFace
1115 (
1116 meshMod,
1117 facei,
1118 masterPointi,
1119 f1,
1120 f1Own,
1121 f1Nei,
1122 patchID
1123 );
1124 }
1125 }
1126 else
1127 {
1128 // f1 is internal face.
1129 if (!modifiedFacei)
1130 {
1131 modFace(meshMod, facei, f1, f1Own, f1Nei, -1);
1132 modifiedFacei = true;
1133 }
1134 else
1135 {
1136 label masterPointi = findPatchFacePoint(f1, -1);
1137
1138 addFace(meshMod, facei, masterPointi, f1, f1Own, f1Nei, -1);
1139 }
1140 }
1141 }
1142
1143 if (f0Own == -1 && f0Nei == -1 && !modifiedFacei)
1144 {
1145 meshMod.setAction(polyRemoveFace(facei));
1146
1147 if (debug & 2)
1148 {
1149 Pout<< "Removed face " << facei << endl;
1150 }
1151 }
1152
1153 faceUptodate[facei] = true;
1154 }
1155
1156
1157 //
1158 // Faces that have not been split but just appended to. Are guaranteed
1159 // to be reachable from an edgeCut.
1160 //
1161
1162 const boolList& edgeIsCut = cuts.edgeIsCut();
1163
1164 forAll(edgeIsCut, edgeI)
1165 {
1166 if (edgeIsCut[edgeI])
1167 {
1168 const labelList& eFaces = mesh().edgeFaces()[edgeI];
1169
1170 forAll(eFaces, i)
1171 {
1172 label facei = eFaces[i];
1173
1174 if (!faceUptodate[facei])
1175 {
1176 // So the face has not been split itself (i.e. its owner
1177 // or neighbour have not been split) so it only
1178 // borders by edge a cell which has been split.
1179
1180 // Get (new or original) owner and neighbour of facei
1181 label own, nei, patchID;
1182 faceCells(cuts, exposedPatchi, facei, own, nei, patchID);
1183
1184
1185 if (own == -1 && nei == -1)
1186 {
1187 meshMod.setAction(polyRemoveFace(facei));
1188
1189 if (debug & 2)
1190 {
1191 Pout<< "Removed face " << facei << endl;
1192 }
1193 }
1194 else
1195 {
1196 // Renumber face to include split edges.
1197 face newFace(addEdgeCutsToFace(facei));
1198
1199 if (debug & 2)
1200 {
1201 Pout<< "Added edge cuts to face " << facei
1202 << " f:" << mesh().faces()[facei]
1203 << " newFace:" << newFace << endl;
1204 }
1205
1206 modFace
1207 (
1208 meshMod,
1209 facei,
1210 newFace,
1211 own,
1212 nei,
1213 patchID
1214 );
1215 }
1216
1217 faceUptodate[facei] = true;
1218 }
1219 }
1220 }
1221 }
1222
1223
1224 //
1225 // Remove any faces on the non-anchor side of a split cell.
1226 // Note: could loop through all cut cells only and check their faces but
1227 // looping over all faces is cleaner and probably faster for dense
1228 // cut patterns.
1229
1230 const faceList& faces = mesh().faces();
1231
1232 forAll(faces, facei)
1233 {
1234 if (!faceUptodate[facei])
1235 {
1236 // Get (new or original) owner and neighbour of facei
1237 label own, nei, patchID;
1238 faceCells(cuts, exposedPatchi, facei, own, nei, patchID);
1239
1240 if (own == -1 && nei == -1)
1241 {
1242 meshMod.setAction(polyRemoveFace(facei));
1243
1244 if (debug & 2)
1245 {
1246 Pout<< "Removed face " << facei << endl;
1247 }
1248 }
1249 else
1250 {
1251 modFace(meshMod, facei, faces[facei], own, nei, patchID);
1252 }
1253
1254 faceUptodate[facei] = true;
1255 }
1256 }
1257
1258 if (debug)
1259 {
1260 Pout<< "meshCutAndRemove:" << nl
1261 << " cells split:" << cuts.nLoops() << nl
1262 << " faces added:" << addedFaces_.size() << nl
1263 << " points added on edges:" << addedPoints_.size() << nl
1264 << endl;
1265 }
1266}
1267
1268
1270{
1271 // Update stored labels for mesh change.
1272 {
1273 Map<label> newAddedFaces(addedFaces_.size());
1274
1275 forAllConstIters(addedFaces_, iter)
1276 {
1277 const label celli = iter.key();
1278 const label addedFacei = iter.val();
1279
1280 const label newCelli = map.reverseCellMap()[celli];
1281 const label newAddedFacei = map.reverseFaceMap()[addedFacei];
1282
1283 if ((newCelli >= 0) && (newAddedFacei >= 0))
1284 {
1285 if
1286 (
1287 (debug & 2)
1288 && (newCelli != celli || newAddedFacei != addedFacei)
1289 )
1290 {
1291 Pout<< "meshCutAndRemove::updateMesh :"
1292 << " updating addedFace for cell " << celli
1293 << " from " << addedFacei
1294 << " to " << newAddedFacei
1295 << endl;
1296 }
1297 newAddedFaces.insert(newCelli, newAddedFacei);
1298 }
1299 }
1300
1301 // Copy
1302 addedFaces_.transfer(newAddedFaces);
1303 }
1304
1305 {
1306 EdgeMap<label> newAddedPoints(addedPoints_.size());
1307
1308 forAllConstIters(addedPoints_, iter)
1309 {
1310 const edge& e = iter.key();
1311 const label addedPointi = iter.val();
1312
1313 const label newStart = map.reversePointMap()[e.start()];
1314 const label newEnd = map.reversePointMap()[e.end()];
1315 const label newAddedPointi = map.reversePointMap()[addedPointi];
1316
1317 if ((newStart >= 0) && (newEnd >= 0) && (newAddedPointi >= 0))
1318 {
1319 edge newE(newStart, newEnd);
1320
1321 if
1322 (
1323 (debug & 2)
1324 && (e != newE || newAddedPointi != addedPointi)
1325 )
1326 {
1327 Pout<< "meshCutAndRemove::updateMesh :"
1328 << " updating addedPoints for edge " << e
1329 << " from " << addedPointi
1330 << " to " << newAddedPointi
1331 << endl;
1332 }
1333
1334 newAddedPoints.insert(newE, newAddedPointi);
1335 }
1336 }
1337
1338 // Copy
1339 addedPoints_.transfer(newAddedPoints);
1340 }
1341}
1342
1343
1344// ************************************************************************* //
Map from edge (expressed as its endpoints) to value. For easier forward declaration it is currently i...
Definition: EdgeMap.H:54
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition: UListI.H:350
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: UList.C:212
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label fcIndex(const label i) const noexcept
Definition: UListI.H:60
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:288
Description of cuts across cells.
Definition: cellCuts.H:113
const labelListList & cellLoops() const
For each cut cell the cut along the circumference.
Definition: cellCuts.H:590
const boolList & pointIsCut() const
Is mesh point cut.
Definition: cellCuts.H:555
const boolList & edgeIsCut() const
Is edge cut.
Definition: cellCuts.H:561
const labelListList & cellAnchorPoints() const
For each cut cell the points on the 'anchor' side of the cell.
Definition: cellCuts.H:602
label nLoops() const
Number of valid cell loops.
Definition: cellCuts.H:596
const scalarField & edgeWeight() const
If edge is cut gives weight (ratio between start() and end())
Definition: cellCuts.H:567
const Map< edge > & faceSplitCut() const
Gives for split face the two cuts that split the face into two.
Definition: cellCuts.H:584
Combines edge or vertex in single label. Used to specify cuts across cell circumference.
Definition: edgeVertex.H:56
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:532
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:469
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:501
Like meshCutter but also removes non-anchor side of cell.
void setRefinement(const label exposedPatchi, const cellCuts &cuts, const labelList &cutPatch, polyTopoChange &meshMod)
Do actual cutting with cut description. Inserts mesh changes.
void updateMesh()
Update for new mesh topology.
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:56
Class containing data for point addition.
Definition: polyAddPoint.H:52
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:498
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1127
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
Class containing data for face removal.
Class containing data for point removal.
Direct mesh changes based on v1.3 polyTopoChange syntax.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const labelListList & cellPoints() const
const labelListList & pointFaces() const
const labelListList & edgeFaces() const
Patchify triangles based on orientation w.r.t other (triangulated or triangulatable) surfaces.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
label nPoints
const labelIOList & zoneID
label findEdge(const edgeList &edges, const labelList &candidates, const label v0, const label v1)
Return edge among candidates that uses the two vertices.
Definition: meshTools.C:359
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition: UListI.H:449
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
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];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } 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};const 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< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
labelList f(nPoints)
labelList pointLabels(nPoints, -1)
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278