cellCuts.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-2017 OpenFOAM Foundation
9 Copyright (C) 2017-2020 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 "cellCuts.H"
30#include "polyMesh.H"
31#include "Time.H"
32#include "ListOps.H"
33#include "cellLooper.H"
34#include "refineCell.H"
35#include "meshTools.H"
36#include "geomCellLooper.H"
37#include "OFstream.H"
38#include "plane.H"
39#include "syncTools.H"
40#include "dummyTransform.H"
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
44namespace Foam
45{
47
48 //- Template specialization for pTraits<edge> so we can use syncTools
49 // functionality
50 template<>
52 {
53 public:
54
55 //- Component type
56 typedef edge cmptType;
57 };
58}
59
60
61// * * * * * * * * * * * * * Private Static Functions * * * * * * * * * * * //
62
63Foam::label Foam::cellCuts::findPartIndex
64(
65 const labelList& elems,
66 const label nElems,
67 const label val
68)
69{
70 for (label i = 0; i < nElems; i++)
71 {
72 if (elems[i] == val)
73 {
74 return i;
75 }
76 }
77 return -1;
78}
79
80
81Foam::boolList Foam::cellCuts::expand
82(
83 const label size,
84 const labelList& labels
85)
86{
87 boolList result(size, false);
88
89 forAll(labels, labelI)
90 {
91 result[labels[labelI]] = true;
92 }
93 return result;
94}
95
96
97Foam::scalarField Foam::cellCuts::expand
98(
99 const label size,
100 const labelList& labels,
101 const scalarField& weights
102)
103{
104 scalarField result(size, -GREAT);
105
106 forAll(labels, labelI)
107 {
108 result[labels[labelI]] = weights[labelI];
109 }
110 return result;
111}
112
113
114Foam::label Foam::cellCuts::firstUnique
115(
116 const labelList& lst,
117 const Map<label>& map
118)
119{
120 forAll(lst, i)
121 {
122 if (!map.found(lst[i]))
123 {
124 return i;
125 }
126 }
127 return -1;
128}
129
130
131// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
132
133void Foam::cellCuts::syncProc()
134{
135 if (!Pstream::parRun())
136 {
137 return;
138 }
139
140 syncTools::syncPointList(mesh(), pointIsCut_, orEqOp<bool>(), false);
141 syncTools::syncEdgeList(mesh(), edgeIsCut_, orEqOp<bool>(), false);
142 syncTools::syncEdgeList(mesh(), edgeWeight_, maxEqOp<scalar>(), -GREAT);
143
144 {
145 const label nBnd = mesh().nBoundaryFaces();
146
147 // Convert faceSplitCut into face-local data: vertex and edge w.r.t.
148 // vertex 0: (since this is same on both sides)
149 //
150 // Sending-side vertex Receiving-side vertex
151 // 0 0
152 // 1 3
153 // 2 2
154 // 3 1
155 //
156 // Sending-side edge Receiving side edge
157 // 0-1 3-0
158 // 1-2 2-3
159 // 2-3 1-2
160 // 3-0 0-1
161 //
162 // Encoding is as index:
163 // 0 : not set
164 // >0 : vertex, vertex is index-1
165 // <0 : edge, edge is -index-1
166
167 edgeList relCuts(nBnd, edge(0, 0));
168
169 const polyBoundaryMesh& pbm = mesh().boundaryMesh();
170
171 forAll(pbm, patchi)
172 {
173 const polyPatch& pp = pbm[patchi];
174
175 if (pp.coupled())
176 {
177 forAll(pp, i)
178 {
179 label facei = pp.start()+i;
180 label bFacei = facei-mesh().nInternalFaces();
181
182 const auto iter = faceSplitCut_.cfind(facei);
183
184 if (iter.found())
185 {
186 const face& f = mesh().faces()[facei];
187 const labelList& fEdges = mesh().faceEdges()[facei];
188 const edge& cuts = iter.val();
189
190 forAll(cuts, i)
191 {
192 if (isEdge(cuts[i]))
193 {
194 label edgei = getEdge(cuts[i]);
195 label index = fEdges.find(edgei);
196 relCuts[bFacei][i] = -index-1;
197 }
198 else
199 {
200 label index = f.find(getVertex(cuts[i]));
201 relCuts[bFacei][i] = index+1;
202 }
203 }
204 }
205 }
206 }
207 }
208
209 // Exchange
211 (
212 mesh(),
213 relCuts,
214 eqOp<edge>(),
215 dummyTransform()
216 );
217
218 // Convert relCuts back into mesh based data
219 forAll(pbm, patchi)
220 {
221 const polyPatch& pp = pbm[patchi];
222
223 if (pp.coupled())
224 {
225 forAll(pp, i)
226 {
227 label facei = pp.start()+i;
228 label bFacei = facei-mesh().nInternalFaces();
229
230 const edge& relCut = relCuts[bFacei];
231 if (relCut != edge(0, 0))
232 {
233 const face& f = mesh().faces()[facei];
234 const labelList& fEdges = mesh().faceEdges()[facei];
235
236 edge absoluteCut(0, 0);
237 forAll(relCut, i)
238 {
239 if (relCut[i] < 0)
240 {
241 label oppFp = -relCut[i]-1;
242 label fp = f.size()-1-oppFp;
243 absoluteCut[i] = edgeToEVert(fEdges[fp]);
244 }
245 else
246 {
247 label oppFp = relCut[i]-1;
248 label fp =
249 (
250 oppFp == 0
251 ? 0
252 : f.size()-oppFp
253 );
254 absoluteCut[i] = vertToEVert(f[fp]);
255 }
256 }
257
258 if
259 (
260 !faceSplitCut_.insert(facei, absoluteCut)
261 && faceSplitCut_[facei] != absoluteCut
262 )
263 {
265 << "Cut " << faceSplitCut_[facei]
266 << " on face " << mesh().faceCentres()[facei]
267 << " of coupled patch " << pp.name()
268 << " is not consistent with coupled cut "
269 << absoluteCut
270 << exit(FatalError);
271 }
272 }
273 }
274 }
275 }
276 }
277}
278
279
280void Foam::cellCuts::writeUncutOBJ
281(
282 const fileName& dir,
283 const label celli
284) const
285{
286 // Cell edges
287 OFstream cutsStream(dir / "cell_" + name(celli) + ".obj");
288
289 Pout<< "Writing cell for time " << mesh().time().timeName()
290 << " to " << cutsStream.name() << nl;
291
293 (
294 cutsStream,
295 mesh().cells(),
296 mesh().faces(),
297 mesh().points(),
298 labelList(1, celli)
299 );
300
301 // Loop cutting cell in two
302 OFstream cutStream(dir / "cellCuts_" + name(celli) + ".obj");
303
304 Pout<< "Writing raw cuts on cell for time " << mesh().time().timeName()
305 << " to " << cutStream.name() << nl;
306
307 const labelList& cPoints = mesh().cellPoints()[celli];
308
309 forAll(cPoints, i)
310 {
311 label pointi = cPoints[i];
312 if (pointIsCut_[pointi])
313 {
314 meshTools::writeOBJ(cutStream, mesh().points()[pointi]);
315 }
316 }
317
318 const pointField& pts = mesh().points();
319
320 const labelList& cEdges = mesh().cellEdges()[celli];
321
322 forAll(cEdges, i)
323 {
324 label edgeI = cEdges[i];
325
326 if (edgeIsCut_[edgeI])
327 {
328 const edge& e = mesh().edges()[edgeI];
329
330 const scalar w = edgeWeight_[edgeI];
331
332 meshTools::writeOBJ(cutStream, w*pts[e[1]] + (1-w)*pts[e[0]]);
333 }
334 }
335}
336
337
338void Foam::cellCuts::writeOBJ
339(
340 const fileName& dir,
341 const label celli,
342 const pointField& loopPoints,
343 const labelList& anchors
344) const
345{
346 // Cell edges
347 OFstream cutsStream(dir / "cell_" + name(celli) + ".obj");
348
349 Pout<< "Writing cell for time " << mesh().time().timeName()
350 << " to " << cutsStream.name() << nl;
351
353 (
354 cutsStream,
355 mesh().cells(),
356 mesh().faces(),
357 mesh().points(),
358 labelList(1, celli)
359 );
360
361
362 // Loop cutting cell in two
363 OFstream loopStream(dir / "cellLoop_" + name(celli) + ".obj");
364
365 Pout<< "Writing loop for time " << mesh().time().timeName()
366 << " to " << loopStream.name() << nl;
367
368 label vertI = 0;
369
370 writeOBJ(loopStream, loopPoints, vertI);
371
372
373 // Anchors for cell
374 OFstream anchorStream(dir / "anchors_" + name(celli) + ".obj");
375
376 Pout<< "Writing anchors for time " << mesh().time().timeName()
377 << " to " << anchorStream.name() << endl;
378
379 forAll(anchors, i)
380 {
381 meshTools::writeOBJ(anchorStream, mesh().points()[anchors[i]]);
382 }
383}
384
385
386Foam::label Foam::cellCuts::edgeEdgeToFace
387(
388 const label celli,
389 const label edgeA,
390 const label edgeB
391) const
392{
393 const labelList& cFaces = mesh().cells()[celli];
394
395 forAll(cFaces, cFacei)
396 {
397 label facei = cFaces[cFacei];
398
399 const labelList& fEdges = mesh().faceEdges()[facei];
400
401 if (fEdges.found(edgeA) && fEdges.found(edgeB))
402 {
403 return facei;
404 }
405 }
406
407 // Coming here means the loop is illegal since the two edges
408 // are not shared by a face. We just mark loop as invalid and continue.
409
411 << "cellCuts : Cannot find face on cell "
412 << celli << " that has both edges " << edgeA << ' ' << edgeB << endl
413 << "faces : " << cFaces << endl
414 << "edgeA : " << mesh().edges()[edgeA] << endl
415 << "edgeB : " << mesh().edges()[edgeB] << endl
416 << "Marking the loop across this cell as invalid" << endl;
417
418 return -1;
419}
420
421
422Foam::label Foam::cellCuts::edgeVertexToFace
423(
424 const label celli,
425 const label edgeI,
426 const label vertI
427) const
428{
429 const labelList& cFaces = mesh().cells()[celli];
430
431 forAll(cFaces, cFacei)
432 {
433 label facei = cFaces[cFacei];
434
435 const face& f = mesh().faces()[facei];
436
437 const labelList& fEdges = mesh().faceEdges()[facei];
438
439 if (fEdges.found(edgeI) && f.found(vertI))
440 {
441 return facei;
442 }
443 }
444
446 << "cellCuts : Cannot find face on cell "
447 << celli << " that has both edge " << edgeI << " and vertex "
448 << vertI << endl
449 << "faces : " << cFaces << endl
450 << "edge : " << mesh().edges()[edgeI] << endl
451 << "Marking the loop across this cell as invalid" << endl;
452
453 return -1;
454}
455
456
457Foam::label Foam::cellCuts::vertexVertexToFace
458(
459 const label celli,
460 const label vertA,
461 const label vertB
462) const
463{
464 const labelList& cFaces = mesh().cells()[celli];
465
466 forAll(cFaces, cFacei)
467 {
468 label facei = cFaces[cFacei];
469
470 const face& f = mesh().faces()[facei];
471
472 if (f.found(vertA) && f.found(vertB))
473 {
474 return facei;
475 }
476 }
477
479 << "cellCuts : Cannot find face on cell "
480 << celli << " that has vertex " << vertA << " and vertex "
481 << vertB << endl
482 << "faces : " << cFaces << endl
483 << "Marking the loop across this cell as invalid" << endl;
484
485 return -1;
486}
487
488
489void Foam::cellCuts::calcFaceCuts() const
490{
491 if (faceCutsPtr_)
492 {
494 << "faceCuts already calculated" << abort(FatalError);
495 }
496
497 const faceList& faces = mesh().faces();
498
499 faceCutsPtr_.reset(new labelListList(mesh().nFaces()));
500 auto& faceCuts = *faceCutsPtr_;
501
502 for (label facei = 0; facei < mesh().nFaces(); facei++)
503 {
504 const face& f = faces[facei];
505
506 // Big enough storage (possibly all points and all edges cut). Shrink
507 // later on.
508 labelList& cuts = faceCuts[facei];
509
510 cuts.setSize(2*f.size());
511
512 label cutI = 0;
513
514 // Do point-edge-point walk over face and collect all cuts.
515 // Problem is that we want to start from one of the endpoints of a
516 // string of connected cuts; we don't want to start somewhere in the
517 // middle.
518
519 // Pass1: find first point cut not preceded by a cut.
520 label startFp = -1;
521
522 forAll(f, fp)
523 {
524 label v0 = f[fp];
525
526 if (pointIsCut_[v0])
527 {
528 label vMin1 = f[f.rcIndex(fp)];
529
530 if
531 (
532 !pointIsCut_[vMin1]
533 && !edgeIsCut_[findEdge(facei, v0, vMin1)]
534 )
535 {
536 cuts[cutI++] = vertToEVert(v0);
537 startFp = f.fcIndex(fp);
538 break;
539 }
540 }
541 }
542
543 // Pass2: first edge cut not preceded by point cut
544 if (startFp == -1)
545 {
546 forAll(f, fp)
547 {
548 label fp1 = f.fcIndex(fp);
549
550 label v0 = f[fp];
551 label v1 = f[fp1];
552
553 label edgeI = findEdge(facei, v0, v1);
554
555 if (edgeIsCut_[edgeI] && !pointIsCut_[v0])
556 {
557 cuts[cutI++] = edgeToEVert(edgeI);
558 startFp = fp1;
559 break;
560 }
561 }
562 }
563
564 // Pass3: nothing found so far. Either face is not cut at all or all
565 // vertices are cut. Start from 0.
566 if (startFp == -1)
567 {
568 startFp = 0;
569 }
570
571 // Store all cuts starting from startFp;
572 label fp = startFp;
573
574 bool allVerticesCut = true;
575
576 forAll(f, i)
577 {
578 label fp1 = f.fcIndex(fp);
579
580 // Get the three items: current vertex, next vertex and edge
581 // inbetween
582 label v0 = f[fp];
583 label v1 = f[fp1];
584 label edgeI = findEdge(facei, v0, v1);
585
586 if (pointIsCut_[v0])
587 {
588 cuts[cutI++] = vertToEVert(v0);
589 }
590 else
591 {
592 // For check if all vertices have been cut (= illegal)
593 allVerticesCut = false;
594 }
595
596 if (edgeIsCut_[edgeI])
597 {
598 cuts[cutI++] = edgeToEVert(edgeI);
599 }
600
601 fp = fp1;
602 }
603
604
605 if (allVerticesCut)
606 {
607 if (verbose_ || debug)
608 {
610 << "Face " << facei << " vertices " << f
611 << " has all its vertices cut. Not cutting face." << endl;
612 }
613
614 cutI = 0;
615 }
616
617 // Remove duplicate starting point
618 if (cutI > 1 && cuts[cutI-1] == cuts[0])
619 {
620 cutI--;
621 }
622 cuts.setSize(cutI);
623 }
624}
625
626
627Foam::label Foam::cellCuts::findEdge
628(
629 const label facei,
630 const label v0,
631 const label v1
632) const
633{
634 const edgeList& edges = mesh().edges();
635
636 const labelList& fEdges = mesh().faceEdges()[facei];
637
638 forAll(fEdges, i)
639 {
640 const edge& e = edges[fEdges[i]];
641
642 if
643 (
644 (e[0] == v0 && e[1] == v1)
645 || (e[0] == v1 && e[1] == v0)
646 )
647 {
648 return fEdges[i];
649 }
650 }
651
652 return -1;
653}
654
655
656Foam::label Foam::cellCuts::loopFace
657(
658 const label celli,
659 const labelList& loop
660) const
661{
662 const cell& cFaces = mesh().cells()[celli];
663
664 forAll(cFaces, cFacei)
665 {
666 label facei = cFaces[cFacei];
667
668 const labelList& fEdges = mesh().faceEdges()[facei];
669 const face& f = mesh().faces()[facei];
670
671 bool allOnFace = true;
672
673 forAll(loop, i)
674 {
675 label cut = loop[i];
676
677 if (isEdge(cut))
678 {
679 if (!fEdges.found(getEdge(cut)))
680 {
681 // Edge not on face. Skip face.
682 allOnFace = false;
683 break;
684 }
685 }
686 else
687 {
688 if (!f.found(getVertex(cut)))
689 {
690 // Vertex not on face. Skip face.
691 allOnFace = false;
692 break;
693 }
694 }
695 }
696
697 if (allOnFace)
698 {
699 // Found face where all elements of loop are on the face.
700 return facei;
701 }
702 }
703 return -1;
704}
705
706
707bool Foam::cellCuts::walkPoint
708(
709 const label celli,
710 const label startCut,
711
712 const label exclude0,
713 const label exclude1,
714
715 const label otherCut,
716
717 label& nVisited,
718 labelList& visited
719) const
720{
721 label vertI = getVertex(otherCut);
722
723 const labelList& pFaces = mesh().pointFaces()[vertI];
724
725 forAll(pFaces, pFacei)
726 {
727 label otherFacei = pFaces[pFacei];
728
729 if
730 (
731 otherFacei != exclude0
732 && otherFacei != exclude1
733 && meshTools::faceOnCell(mesh(), celli, otherFacei)
734 )
735 {
736 label oldNVisited = nVisited;
737
738 bool foundLoop =
739 walkCell
740 (
741 celli,
742 startCut,
743 otherFacei,
744 otherCut,
745 nVisited,
746 visited
747 );
748
749 if (foundLoop)
750 {
751 return true;
752 }
753
754 // No success. Restore state and continue
755 nVisited = oldNVisited;
756 }
757 }
758 return false;
759}
760
761
762bool Foam::cellCuts::crossEdge
763(
764 const label celli,
765 const label startCut,
766 const label facei,
767 const label otherCut,
768
769 label& nVisited,
770 labelList& visited
771) const
772{
773 // Cross edge to other face
774 label edgeI = getEdge(otherCut);
775
776 label otherFacei = meshTools::otherFace(mesh(), celli, facei, edgeI);
777
778 // Store old state
779 label oldNVisited = nVisited;
780
781 // Recurse to otherCut
782 bool foundLoop =
783 walkCell
784 (
785 celli,
786 startCut,
787 otherFacei,
788 otherCut,
789 nVisited,
790 visited
791 );
792
793 if (foundLoop)
794 {
795 return true;
796 }
797 else
798 {
799 // No success. Restore state (i.e. backtrack)
800 nVisited = oldNVisited;
801
802 return false;
803 }
804}
805
806
807bool Foam::cellCuts::addCut
808(
809 const label celli,
810 const label cut,
811 label& nVisited,
812 labelList& visited
813) const
814{
815 // Check for duplicate cuts.
816 if (findPartIndex(visited, nVisited, cut) != -1)
817 {
818 // Truncate (copy of) visited for ease of printing.
819 labelList truncVisited(visited);
820 truncVisited.setSize(nVisited);
821
822 if (verbose_ || debug)
823 {
824 Pout<< "For cell " << celli << " : trying to add duplicate cut "
825 << cut;
826 labelList cuts(1, cut);
827 writeCuts(Pout, cuts, loopWeights(cuts));
828
829 Pout<< " to path:";
830 writeCuts(Pout, truncVisited, loopWeights(truncVisited));
831 Pout<< endl;
832 }
833
834 return false;
835 }
836 else
837 {
838 visited[nVisited++] = cut;
839
840 return true;
841 }
842}
843
844
845bool Foam::cellCuts::walkFace
846(
847 const label celli,
848 const label startCut,
849 const label facei,
850 const label cut,
851
852 label& lastCut,
853 label& beforeLastCut,
854 label& nVisited,
855 labelList& visited
856) const
857{
858 const labelList& fCuts = faceCuts()[facei];
859
860 if (fCuts.size() < 2)
861 {
862 return false;
863 }
864
865 // Easy case : two cuts.
866 if (fCuts.size() == 2)
867 {
868 if (fCuts[0] == cut)
869 {
870 if (!addCut(celli, cut, nVisited, visited))
871 {
872 return false;
873 }
874
875 beforeLastCut = cut;
876 lastCut = fCuts[1];
877
878 return true;
879 }
880 else
881 {
882 if (!addCut(celli, cut, nVisited, visited))
883 {
884 return false;
885 }
886
887 beforeLastCut = cut;
888 lastCut = fCuts[0];
889
890 return true;
891 }
892 }
893
894 // Harder case: more than 2 cuts on face.
895 // Should be start or end of string of cuts. Store all but last cut.
896 if (fCuts[0] == cut)
897 {
898 // Walk forward
899 for (label i = 0; i < fCuts.size()-1; i++)
900 {
901 if (!addCut(celli, fCuts[i], nVisited, visited))
902 {
903 return false;
904 }
905 }
906 beforeLastCut = fCuts[fCuts.size()-2];
907 lastCut = fCuts[fCuts.size()-1];
908 }
909 else if (fCuts[fCuts.size()-1] == cut)
910 {
911 for (label i = fCuts.size()-1; i >= 1; --i)
912 {
913 if (!addCut(celli, fCuts[i], nVisited, visited))
914 {
915 return false;
916 }
917 }
918 beforeLastCut = fCuts[1];
919 lastCut = fCuts[0];
920 }
921 else
922 {
923 if (verbose_ || debug)
924 {
926 << "In middle of cut. cell:" << celli << " face:" << facei
927 << " cuts:" << fCuts << " current cut:" << cut << endl;
928 }
929
930 return false;
931 }
932
933 return true;
934}
935
936
937
938bool Foam::cellCuts::walkCell
939(
940 const label celli,
941 const label startCut,
942 const label facei,
943 const label cut,
944
945 label& nVisited,
946 labelList& visited
947) const
948{
949 // Walk across face, storing cuts. Return last two cuts visited.
950 label lastCut = -1;
951 label beforeLastCut = -1;
952
953
954 if (debug & 2)
955 {
956 Pout<< "For cell:" << celli << " walked across face " << facei
957 << " from cut ";
958 labelList cuts(1, cut);
959 writeCuts(Pout, cuts, loopWeights(cuts));
960 Pout<< endl;
961 }
962
963 bool validWalk = walkFace
964 (
965 celli,
966 startCut,
967 facei,
968 cut,
969
970 lastCut,
971 beforeLastCut,
972 nVisited,
973 visited
974 );
975
976 if (!validWalk)
977 {
978 return false;
979 }
980
981 if (debug & 2)
982 {
983 Pout<< " to last cut ";
984 labelList cuts(1, lastCut);
985 writeCuts(Pout, cuts, loopWeights(cuts));
986 Pout<< endl;
987 }
988
989 // Check if starting point reached.
990 if (lastCut == startCut)
991 {
992 if (nVisited >= 3)
993 {
994 if (debug & 2)
995 {
996 // Truncate (copy of) visited for ease of printing.
997 labelList truncVisited(visited);
998 truncVisited.setSize(nVisited);
999
1000 Pout<< "For cell " << celli << " : found closed path:";
1001 writeCuts(Pout, truncVisited, loopWeights(truncVisited));
1002 Pout<< " closed by " << lastCut << endl;
1003 }
1004
1005 return true;
1006 }
1007 else
1008 {
1009 // Loop but too short. Probably folded back on itself.
1010 return false;
1011 }
1012 }
1013
1014
1015 // Check last cut and one before that to determine type.
1016
1017 // From beforeLastCut to lastCut:
1018 // - from edge to edge
1019 // (- always not along existing edge)
1020 // - from edge to vertex
1021 // - not along existing edge
1022 // - along existing edge: not allowed?
1023 // - from vertex to edge
1024 // - not along existing edge
1025 // - along existing edge. Not allowed. See above.
1026 // - from vertex to vertex
1027 // - not along existing edge
1028 // - along existing edge
1029
1030 if (isEdge(beforeLastCut))
1031 {
1032 if (isEdge(lastCut))
1033 {
1034 // beforeLastCut=edge, lastCut=edge.
1035
1036 // Cross lastCut (=edge) into face not facei
1037 return crossEdge
1038 (
1039 celli,
1040 startCut,
1041 facei,
1042 lastCut,
1043 nVisited,
1044 visited
1045 );
1046 }
1047 else
1048 {
1049 // beforeLastCut=edge, lastCut=vertex.
1050
1051 // Go from lastCut to all connected faces (except facei)
1052 return walkPoint
1053 (
1054 celli,
1055 startCut,
1056 facei, // exclude0: don't cross back on current face
1057 -1, // exclude1
1058 lastCut,
1059 nVisited,
1060 visited
1061 );
1062 }
1063 }
1064 else
1065 {
1066 if (isEdge(lastCut))
1067 {
1068 // beforeLastCut=vertex, lastCut=edge.
1069 return crossEdge
1070 (
1071 celli,
1072 startCut,
1073 facei,
1074 lastCut,
1075 nVisited,
1076 visited
1077 );
1078 }
1079 else
1080 {
1081 // beforeLastCut=vertex, lastCut=vertex. Check if along existing
1082 // edge.
1083 label edgeI =
1084 findEdge
1085 (
1086 facei,
1087 getVertex(beforeLastCut),
1088 getVertex(lastCut)
1089 );
1090
1091 if (edgeI != -1)
1092 {
1093 // Cut along existing edge. So is in fact on two faces.
1094 // Get faces on both sides of the edge to make
1095 // sure we do not fold back on to those.
1096
1097 label f0, f1;
1098 meshTools::getEdgeFaces(mesh(), celli, edgeI, f0, f1);
1099
1100 return walkPoint
1101 (
1102 celli,
1103 startCut,
1104 f0,
1105 f1,
1106 lastCut,
1107 nVisited,
1108 visited
1109 );
1110
1111 }
1112 else
1113 {
1114 // Cross cut across face.
1115 return walkPoint
1116 (
1117 celli,
1118 startCut,
1119 facei, // face to exclude
1120 -1, // face to exclude
1121 lastCut,
1122 nVisited,
1123 visited
1124 );
1125 }
1126 }
1127 }
1128}
1129
1130
1131void Foam::cellCuts::calcCellLoops(const labelList& cutCells)
1132{
1133 // Determine for every cut cell the loop (= face) it is cut by. Done by
1134 // starting from a cut edge or cut vertex and walking across faces, from
1135 // cut to cut, until starting cut hit.
1136 // If multiple loops are possible across a cell circumference takes the
1137 // first one found.
1138
1139 // Calculate cuts per face.
1140 const labelListList& allFaceCuts = faceCuts();
1141
1142 // Per cell the number of faces with valid cuts. Is used as quick
1143 // rejection to see if cell can be cut.
1144 labelList nCutFaces(mesh().nCells(), Zero);
1145
1146 forAll(allFaceCuts, facei)
1147 {
1148 const labelList& fCuts = allFaceCuts[facei];
1149
1150 if (fCuts.size() == mesh().faces()[facei].size())
1151 {
1152 // Too many cuts on face. WalkCell would get very upset so disable.
1153 nCutFaces[mesh().faceOwner()[facei]] = labelMin;
1154
1155 if (mesh().isInternalFace(facei))
1156 {
1157 nCutFaces[mesh().faceNeighbour()[facei]] = labelMin;
1158 }
1159 }
1160 else if (fCuts.size() >= 2)
1161 {
1162 // Could be valid cut. Update count for owner and neighbour.
1163 nCutFaces[mesh().faceOwner()[facei]]++;
1164
1165 if (mesh().isInternalFace(facei))
1166 {
1167 nCutFaces[mesh().faceNeighbour()[facei]]++;
1168 }
1169 }
1170 }
1171
1172
1173 // Stack of visited cuts (nVisited used as stack pointer)
1174 // Size big enough.
1175 labelList visited(mesh().nPoints());
1176
1177 forAll(cutCells, i)
1178 {
1179 label celli = cutCells[i];
1180
1181 bool validLoop = false;
1182
1183 // Quick rejection: has enough faces that are cut?
1184 if (nCutFaces[celli] >= 1)
1185 {
1186 const labelList& cFaces = mesh().cells()[celli];
1187
1188 if (debug & 2)
1189 {
1190 Pout<< "cell:" << celli << " cut faces:" << endl;
1191 forAll(cFaces, i)
1192 {
1193 label facei = cFaces[i];
1194 const labelList& fCuts = allFaceCuts[facei];
1195
1196 Pout<< " face:" << facei << " cuts:";
1197 writeCuts(Pout, fCuts, loopWeights(fCuts));
1198 Pout<< endl;
1199 }
1200 }
1201
1202 label nVisited = 0;
1203
1204 // Determine the first cut face to start walking from.
1205 forAll(cFaces, cFacei)
1206 {
1207 label facei = cFaces[cFacei];
1208
1209 const labelList& fCuts = allFaceCuts[facei];
1210
1211 // Take first or last cut of multiple on face.
1212 // Note that in calcFaceCuts
1213 // we have already made sure this is the start or end of a
1214 // string of cuts.
1215 if (fCuts.size() >= 2)
1216 {
1217 // Try walking from start of fCuts.
1218 nVisited = 0;
1219
1220 if (debug & 2)
1221 {
1222 Pout<< "cell:" << celli
1223 << " start walk at face:" << facei
1224 << " cut:";
1225 labelList cuts(1, fCuts[0]);
1226 writeCuts(Pout, cuts, loopWeights(cuts));
1227 Pout<< endl;
1228 }
1229
1230 validLoop =
1231 walkCell
1232 (
1233 celli,
1234 fCuts[0],
1235 facei,
1236 fCuts[0],
1237
1238 nVisited,
1239 visited
1240 );
1241
1242 if (validLoop)
1243 {
1244 break;
1245 }
1246
1247 // No need to try and walk from end of cuts since
1248 // all paths already tried by walkCell.
1249 }
1250 }
1251
1252 if (validLoop)
1253 {
1254 // Copy nVisited elements out of visited (since visited is
1255 // never truncated for efficiency reasons)
1256
1257 labelList& loop = cellLoops_[celli];
1258
1259 loop.setSize(nVisited);
1260
1261 for (label i = 0; i < nVisited; i++)
1262 {
1263 loop[i] = visited[i];
1264 }
1265 }
1266 else
1267 {
1268 // Invalid loop. Leave cellLoops_[celli] zero size which
1269 // flags this.
1270 if (verbose_ || debug)
1271 {
1272 Pout<< "calcCellLoops(const labelList&) :"
1273 << " did not find valid"
1274 << " loop for cell " << celli << endl;
1275 // Dump cell and cuts on cell.
1276 writeUncutOBJ(".", celli);
1277 }
1278 cellLoops_[celli].clear();
1279 }
1280 }
1281 else
1282 {
1283 //Pout<< "calcCellLoops(const labelList&) : did not find valid"
1284 // << " loop for cell " << celli << " since not enough cut faces"
1285 // << endl;
1286 cellLoops_[celli].clear();
1287 }
1288 }
1289}
1290
1291
1292void Foam::cellCuts::walkEdges
1293(
1294 const label celli,
1295 const label pointi,
1296 const label status,
1297
1298 Map<label>& edgeStatus,
1299 Map<label>& pointStatus
1300) const
1301{
1302 if (pointStatus.insert(pointi, status))
1303 {
1304 // First visit to pointi
1305
1306 const labelList& pEdges = mesh().pointEdges()[pointi];
1307
1308 forAll(pEdges, pEdgeI)
1309 {
1310 label edgeI = pEdges[pEdgeI];
1311
1312 if
1313 (
1314 meshTools::edgeOnCell(mesh(), celli, edgeI)
1315 && edgeStatus.insert(edgeI, status)
1316 )
1317 {
1318 // First visit to edgeI so recurse.
1319
1320 label v2 = mesh().edges()[edgeI].otherVertex(pointi);
1321
1322 walkEdges(celli, v2, status, edgeStatus, pointStatus);
1323 }
1324 }
1325 }
1326}
1327
1328
1330(
1331 const labelList& cellPoints,
1332 const labelList& anchorPoints,
1333 const labelList& loop
1334) const
1335{
1336 labelList newElems(cellPoints.size());
1337 label newElemI = 0;
1338
1339 forAll(cellPoints, i)
1340 {
1341 label pointi = cellPoints[i];
1342
1343 if
1344 (
1345 !anchorPoints.found(pointi)
1346 && !loop.found(vertToEVert(pointi))
1347 )
1348 {
1349 newElems[newElemI++] = pointi;
1350 }
1351 }
1352
1353 newElems.setSize(newElemI);
1354
1355 return newElems;
1356}
1357
1358
1359bool Foam::cellCuts::loopAnchorConsistent
1360(
1361 const label celli,
1362 const pointField& loopPts,
1363 const labelList& anchorPoints
1364) const
1365{
1366 // Create identity face for ease of calculation of normal etc.
1367 face f(identity(loopPts.size()));
1368
1369 const vector areaNorm = f.areaNormal(loopPts);
1370 const point ctr = f.centre(loopPts);
1371
1372
1373 // Get average position of anchor points.
1374 vector avg(Zero);
1375
1376 for (const label pointi : anchorPoints)
1377 {
1378 avg += mesh().points()[pointi];
1379 }
1380 avg /= anchorPoints.size();
1381
1382
1383 if (((avg - ctr) & areaNorm) > 0)
1384 {
1385 return true;
1386 }
1387
1388 return false;
1389}
1390
1391
1392bool Foam::cellCuts::calcAnchors
1393(
1394 const label celli,
1395 const labelList& loop,
1396 const pointField& loopPts,
1397
1398 labelList& anchorPoints
1399) const
1400{
1401 const edgeList& edges = mesh().edges();
1402
1403 const labelList& cPoints = mesh().cellPoints()[celli];
1404 const labelList& cEdges = mesh().cellEdges()[celli];
1405 const cell& cFaces = mesh().cells()[celli];
1406
1407 // Points on loop
1408
1409 // Status of point:
1410 // 0 - on loop
1411 // 1 - point set 1
1412 // 2 - point set 2
1413 Map<label> pointStatus(2*cPoints.size());
1414 Map<label> edgeStatus(2*cEdges.size());
1415
1416 // Mark loop vertices
1417 forAll(loop, i)
1418 {
1419 label cut = loop[i];
1420
1421 if (isEdge(cut))
1422 {
1423 edgeStatus.insert(getEdge(cut), 0);
1424 }
1425 else
1426 {
1427 pointStatus.insert(getVertex(cut), 0);
1428 }
1429 }
1430 // Since edges between two cut vertices have not been marked, mark them
1431 // explicitly
1432 forAll(cEdges, i)
1433 {
1434 label edgeI = cEdges[i];
1435 const edge& e = edges[edgeI];
1436
1437 if (pointStatus.found(e[0]) && pointStatus.found(e[1]))
1438 {
1439 edgeStatus.insert(edgeI, 0);
1440 }
1441 }
1442
1443
1444 // Find uncut starting vertex
1445 label uncutIndex = firstUnique(cPoints, pointStatus);
1446
1447 if (uncutIndex == -1)
1448 {
1449 if (verbose_ || debug)
1450 {
1452 << "Invalid loop " << loop << " for cell " << celli << endl
1453 << "Can not find point on cell which is not cut by loop."
1454 << endl;
1455
1456 writeOBJ(".", celli, loopPts, labelList(0));
1457 }
1458
1459 return false;
1460 }
1461
1462 // Walk unset vertices and edges and mark with 1 in pointStatus, edgeStatus
1463 walkEdges(celli, cPoints[uncutIndex], 1, edgeStatus, pointStatus);
1464
1465 // Find new uncut starting vertex
1466 uncutIndex = firstUnique(cPoints, pointStatus);
1467
1468 if (uncutIndex == -1)
1469 {
1470 // All vertices either in loop or in anchor. So split is along single
1471 // face.
1472 if (verbose_ || debug)
1473 {
1475 << "Invalid loop " << loop << " for cell " << celli << endl
1476 << "All vertices of cell are either in loop or in anchor set"
1477 << endl;
1478
1479 writeOBJ(".", celli, loopPts, labelList(0));
1480 }
1481
1482 return false;
1483 }
1484
1485 // Walk unset vertices and edges and mark with 2. These are the
1486 // pointset 2.
1487 walkEdges(celli, cPoints[uncutIndex], 2, edgeStatus, pointStatus);
1488
1489 // Collect both sets in lists.
1490 DynamicList<label> connectedPoints(cPoints.size());
1491 DynamicList<label> otherPoints(cPoints.size());
1492
1493 forAllConstIters(pointStatus, iter)
1494 {
1495 if (iter.val() == 1)
1496 {
1497 connectedPoints.append(iter.key());
1498 }
1499 else if (iter.val() == 2)
1500 {
1501 otherPoints.append(iter.key());
1502 }
1503 }
1504 connectedPoints.shrink();
1505 otherPoints.shrink();
1506
1507 // Check that all points have been used.
1508 uncutIndex = firstUnique(cPoints, pointStatus);
1509
1510 if (uncutIndex != -1)
1511 {
1512 if (verbose_ || debug)
1513 {
1515 << "Invalid loop " << loop << " for cell " << celli
1516 << " since it splits the cell into more than two cells" << endl;
1517
1518 writeOBJ(".", celli, loopPts, connectedPoints);
1519 }
1520
1521 return false;
1522 }
1523
1524
1525 // Check that both parts (connectedPoints, otherPoints) have enough faces.
1526 labelHashSet connectedFaces(2*cFaces.size());
1527 labelHashSet otherFaces(2*cFaces.size());
1528
1529 forAllConstIters(pointStatus, iter)
1530 {
1531 const label pointi = iter.key();
1532
1533 const labelList& pFaces = mesh().pointFaces()[pointi];
1534
1535 if (iter.val() == 1)
1536 {
1537 forAll(pFaces, pFacei)
1538 {
1539 if (meshTools::faceOnCell(mesh(), celli, pFaces[pFacei]))
1540 {
1541 connectedFaces.insert(pFaces[pFacei]);
1542 }
1543 }
1544 }
1545 else if (iter.val() == 2)
1546 {
1547 forAll(pFaces, pFacei)
1548 {
1549 if (meshTools::faceOnCell(mesh(), celli, pFaces[pFacei]))
1550 {
1551 otherFaces.insert(pFaces[pFacei]);
1552 }
1553 }
1554 }
1555 }
1556
1557 if (connectedFaces.size() < 3)
1558 {
1559 if (verbose_ || debug)
1560 {
1562 << "Invalid loop " << loop << " for cell " << celli
1563 << " since would have too few faces on one side." << nl
1564 << "All faces:" << cFaces << endl;
1565
1566 writeOBJ(".", celli, loopPts, connectedPoints);
1567 }
1568
1569 return false;
1570 }
1571
1572 if (otherFaces.size() < 3)
1573 {
1574 if (verbose_ || debug)
1575 {
1577 << "Invalid loop " << loop << " for cell " << celli
1578 << " since would have too few faces on one side." << nl
1579 << "All faces:" << cFaces << endl;
1580
1581 writeOBJ(".", celli, loopPts, otherPoints);
1582 }
1583
1584 return false;
1585 }
1586
1587
1588 // Check that faces are split into two regions and not more.
1589 // When walking across the face the only transition of pointStatus is
1590 // from set1 to loop to set2 or back. Not allowed is from set1 to loop to
1591 // set1.
1592 {
1593 forAll(cFaces, i)
1594 {
1595 label facei = cFaces[i];
1596
1597 const face& f = mesh().faces()[facei];
1598
1599 bool hasSet1 = false;
1600 bool hasSet2 = false;
1601
1602 label prevStat = pointStatus[f[0]];
1603
1604 forAll(f, fp)
1605 {
1606 label v0 = f[fp];
1607 label pStat = pointStatus[v0];
1608
1609 if (pStat == prevStat)
1610 {
1611 }
1612 else if (pStat == 0)
1613 {
1614 // Loop.
1615 }
1616 else if (pStat == 1)
1617 {
1618 if (hasSet1)
1619 {
1620 // Second occurrence of set1.
1621 if (verbose_ || debug)
1622 {
1624 << "Invalid loop " << loop << " for cell "
1625 << celli
1626 << " since face " << f << " would be split into"
1627 << " more than two faces" << endl;
1628
1629 writeOBJ(".", celli, loopPts, otherPoints);
1630 }
1631 return false;
1632 }
1633
1634 hasSet1 = true;
1635 }
1636 else if (pStat == 2)
1637 {
1638 if (hasSet2)
1639 {
1640 // Second occurrence of set1.
1641 if (verbose_ || debug)
1642 {
1644 << "Invalid loop " << loop << " for cell "
1645 << celli
1646 << " since face " << f << " would be split into"
1647 << " more than two faces" << endl;
1648
1649 writeOBJ(".", celli, loopPts, otherPoints);
1650 }
1651 return false;
1652 }
1653
1654 hasSet2 = true;
1655 }
1656 else
1657 {
1659 << abort(FatalError);
1660 }
1661
1662 prevStat = pStat;
1663
1664
1665 label v1 = f.nextLabel(fp);
1666 label edgeI = findEdge(facei, v0, v1);
1667
1668 label eStat = edgeStatus[edgeI];
1669
1670 if (eStat == prevStat)
1671 {
1672 }
1673 else if (eStat == 0)
1674 {
1675 // Loop.
1676 }
1677 else if (eStat == 1)
1678 {
1679 if (hasSet1)
1680 {
1681 // Second occurrence of set1.
1682 if (verbose_ || debug)
1683 {
1685 << "Invalid loop " << loop << " for cell "
1686 << celli
1687 << " since face " << f << " would be split into"
1688 << " more than two faces" << endl;
1689
1690 writeOBJ(".", celli, loopPts, otherPoints);
1691 }
1692
1693 return false;
1694 }
1695
1696 hasSet1 = true;
1697 }
1698 else if (eStat == 2)
1699 {
1700 if (hasSet2)
1701 {
1702 // Second occurrence of set1.
1703 if (verbose_ || debug)
1704 {
1706 << "Invalid loop " << loop << " for cell "
1707 << celli
1708 << " since face " << f << " would be split into"
1709 << " more than two faces" << endl;
1710
1711 writeOBJ(".", celli, loopPts, otherPoints);
1712 }
1713 return false;
1714 }
1715
1716 hasSet2 = true;
1717 }
1718 prevStat = eStat;
1719 }
1720 }
1721 }
1722
1723
1724
1725
1726 // Check which one of point sets to use.
1727 bool loopOk = loopAnchorConsistent(celli, loopPts, connectedPoints);
1728
1729 //if (debug)
1730 {
1731 // Additional check: are non-anchor points on other side?
1732 bool otherLoopOk = loopAnchorConsistent(celli, loopPts, otherPoints);
1733
1734 if (loopOk == otherLoopOk)
1735 {
1736 // Both sets of points are supposedly on the same side as the
1737 // loop normal. Oops.
1738
1740 << " For cell:" << celli
1741 << " achorpoints and nonanchorpoints are geometrically"
1742 << " on same side!" << endl
1743 << "cellPoints:" << cPoints << endl
1744 << "loop:" << loop << endl
1745 << "anchors:" << connectedPoints << endl
1746 << "otherPoints:" << otherPoints << endl;
1747
1748 writeOBJ(".", celli, loopPts, connectedPoints);
1749 }
1750 }
1751
1752 if (loopOk)
1753 {
1754 // connectedPoints on 'outside' of loop so these are anchor points
1755 anchorPoints.transfer(connectedPoints);
1756 connectedPoints.clear();
1757 }
1758 else
1759 {
1760 anchorPoints.transfer(otherPoints);
1761 otherPoints.clear();
1762 }
1763 return true;
1764}
1765
1766
1767Foam::pointField Foam::cellCuts::loopPoints
1768(
1769 const labelList& loop,
1770 const scalarField& loopWeights
1771) const
1772{
1773 pointField loopPts(loop.size());
1774
1775 forAll(loop, fp)
1776 {
1777 loopPts[fp] = coord(loop[fp], loopWeights[fp]);
1778 }
1779 return loopPts;
1780}
1781
1782
1783Foam::scalarField Foam::cellCuts::loopWeights(const labelList& loop) const
1784{
1785 scalarField weights(loop.size());
1786
1787 forAll(loop, fp)
1788 {
1789 label cut = loop[fp];
1790
1791 if (isEdge(cut))
1792 {
1793 label edgeI = getEdge(cut);
1794
1795 weights[fp] = edgeWeight_[edgeI];
1796 }
1797 else
1798 {
1799 weights[fp] = -GREAT;
1800 }
1801 }
1802 return weights;
1803}
1804
1805
1806bool Foam::cellCuts::validEdgeLoop
1807(
1808 const labelList& loop,
1809 const scalarField& loopWeights
1810) const
1811{
1812 forAll(loop, fp)
1813 {
1814 label cut = loop[fp];
1815
1816 if (isEdge(cut))
1817 {
1818 label edgeI = getEdge(cut);
1819
1820 // Check: cut compatible only if can be snapped to existing one.
1821 if (edgeIsCut_[edgeI])
1822 {
1823 scalar edgeLen = mesh().edges()[edgeI].mag(mesh().points());
1824
1825 if
1826 (
1827 mag(loopWeights[fp] - edgeWeight_[edgeI])
1828 > geomCellLooper::snapTol()*edgeLen
1829 )
1830 {
1831 // Positions differ too much->would create two cuts.
1832 return false;
1833 }
1834 }
1835 }
1836 }
1837 return true;
1838}
1839
1840
1841Foam::label Foam::cellCuts::countFaceCuts
1842(
1843 const label facei,
1844 const labelList& loop
1845) const
1846{
1847 // Includes cuts through vertices and through edges.
1848 // Assumes that if edge is cut both in edgeIsCut and in loop that the
1849 // position of the cut is the same.
1850
1851 label nCuts = 0;
1852
1853 // Count cut vertices
1854 const face& f = mesh().faces()[facei];
1855
1856 forAll(f, fp)
1857 {
1858 label vertI = f[fp];
1859
1860 // Vertex already cut or mentioned in current loop.
1861 if
1862 (
1863 pointIsCut_[vertI]
1864 || loop.found(vertToEVert(vertI))
1865 )
1866 {
1867 nCuts++;
1868 }
1869 }
1870
1871 // Count cut edges.
1872 const labelList& fEdges = mesh().faceEdges()[facei];
1873
1874 forAll(fEdges, fEdgeI)
1875 {
1876 label edgeI = fEdges[fEdgeI];
1877
1878 // Edge already cut or mentioned in current loop.
1879 if
1880 (
1881 edgeIsCut_[edgeI]
1882 || loop.found(edgeToEVert(edgeI))
1883 )
1884 {
1885 nCuts++;
1886 }
1887 }
1888
1889 return nCuts;
1890}
1891
1892
1893bool Foam::cellCuts::conservativeValidLoop
1894(
1895 const label celli,
1896 const labelList& loop
1897) const
1898{
1899
1900 if (loop.size() < 2)
1901 {
1902 return false;
1903 }
1904
1905 forAll(loop, cutI)
1906 {
1907 if (isEdge(loop[cutI]))
1908 {
1909 label edgeI = getEdge(loop[cutI]);
1910
1911 if (edgeIsCut_[edgeI])
1912 {
1913 // edge compatibility already checked.
1914 }
1915 else
1916 {
1917 // Quick rejection: vertices of edge itself cannot be cut.
1918 const edge& e = mesh().edges()[edgeI];
1919
1920 if (pointIsCut_[e.start()] || pointIsCut_[e.end()])
1921 {
1922 return false;
1923 }
1924
1925
1926 // Check faces using this edge
1927 const labelList& eFaces = mesh().edgeFaces()[edgeI];
1928
1929 forAll(eFaces, eFacei)
1930 {
1931 label nCuts = countFaceCuts(eFaces[eFacei], loop);
1932
1933 if (nCuts > 2)
1934 {
1935 return false;
1936 }
1937 }
1938 }
1939 }
1940 else
1941 {
1942 // Vertex cut
1943
1944 label vertI = getVertex(loop[cutI]);
1945
1946 if (!pointIsCut_[vertI])
1947 {
1948 // New cut through vertex.
1949
1950 // Check edges using vertex.
1951 const labelList& pEdges = mesh().pointEdges()[vertI];
1952
1953 forAll(pEdges, pEdgeI)
1954 {
1955 label edgeI = pEdges[pEdgeI];
1956
1957 if (edgeIsCut_[edgeI])
1958 {
1959 return false;
1960 }
1961 }
1962
1963 // Check faces using vertex.
1964 const labelList& pFaces = mesh().pointFaces()[vertI];
1965
1966 forAll(pFaces, pFacei)
1967 {
1968 label nCuts = countFaceCuts(pFaces[pFacei], loop);
1969
1970 if (nCuts > 2)
1971 {
1972 return false;
1973 }
1974 }
1975 }
1976 }
1977 }
1978
1979 // All ok.
1980 return true;
1981}
1982
1983
1984bool Foam::cellCuts::validLoop
1985(
1986 const label celli,
1987 const labelList& loop,
1988 const scalarField& loopWeights,
1989
1990 Map<edge>& newFaceSplitCut,
1991 labelList& anchorPoints
1992) const
1993{
1994 // Determine compatibility of loop with existing cut pattern. Does not use
1995 // derived cut-addressing (faceCuts), only pointIsCut, edgeIsCut.
1996 // Adds any cross-cuts found to newFaceSplitCut and sets cell points on
1997 // one side of the loop in anchorPoints.
1998
1999 if (loop.size() < 2)
2000 {
2001 return false;
2002 }
2003
2004 if (debug & 4)
2005 {
2006 // Allow as fallback the 'old' loop checking where only a single
2007 // cut per face is allowed.
2008 if (!conservativeValidLoop(celli, loop))
2009 {
2010 Info << "Invalid conservative loop: " << loop << endl;
2011 return false;
2012 }
2013 }
2014
2015 forAll(loop, fp)
2016 {
2017 label cut = loop[fp];
2018 label nextCut = loop[(fp+1) % loop.size()];
2019
2020 // Label (if any) of face cut (so cut not along existing edge)
2021 label meshFacei = -1;
2022
2023 if (isEdge(cut))
2024 {
2025 label edgeI = getEdge(cut);
2026
2027 // Look one cut ahead to find if it is along existing edge.
2028
2029 if (isEdge(nextCut))
2030 {
2031 // From edge to edge -> cross cut
2032 label nextEdgeI = getEdge(nextCut);
2033
2034 // Find face and mark as to be split.
2035 meshFacei = edgeEdgeToFace(celli, edgeI, nextEdgeI);
2036
2037 if (meshFacei == -1)
2038 {
2039 // Can't find face using both cut edges.
2040 return false;
2041 }
2042 }
2043 else
2044 {
2045 // From edge to vertex -> cross cut only if no existing edge.
2046
2047 label nextVertI = getVertex(nextCut);
2048 const edge& e = mesh().edges()[edgeI];
2049
2050 if (e.start() != nextVertI && e.end() != nextVertI)
2051 {
2052 // New edge. Find face and mark as to be split.
2053 meshFacei = edgeVertexToFace(celli, edgeI, nextVertI);
2054
2055 if (meshFacei == -1)
2056 {
2057 // Can't find face. Illegal.
2058 return false;
2059 }
2060 }
2061 }
2062 }
2063 else
2064 {
2065 // Cut is vertex.
2066 label vertI = getVertex(cut);
2067
2068 if (isEdge(nextCut))
2069 {
2070 // From vertex to edge -> cross cut only if no existing edge.
2071 label nextEdgeI = getEdge(nextCut);
2072
2073 const edge& nextE = mesh().edges()[nextEdgeI];
2074
2075 if (nextE.start() != vertI && nextE.end() != vertI)
2076 {
2077 // Should be cross cut. Find face.
2078 meshFacei = edgeVertexToFace(celli, nextEdgeI, vertI);
2079
2080 if (meshFacei == -1)
2081 {
2082 return false;
2083 }
2084 }
2085 }
2086 else
2087 {
2088 // From vertex to vertex -> cross cut only if no existing edge.
2089 label nextVertI = getVertex(nextCut);
2090
2091 if (meshTools::findEdge(mesh(), vertI, nextVertI) == -1)
2092 {
2093 // New edge. Find face.
2094 meshFacei = vertexVertexToFace(celli, vertI, nextVertI);
2095
2096 if (meshFacei == -1)
2097 {
2098 return false;
2099 }
2100 }
2101 }
2102 }
2103
2104 if (meshFacei != -1)
2105 {
2106 // meshFacei is cut across along cut-nextCut (so not along existing
2107 // edge). Check if this is compatible with existing pattern.
2108 edge cutEdge(cut, nextCut);
2109
2110 const auto iter = faceSplitCut_.cfind(meshFacei);
2111
2112 if (!iter.found())
2113 {
2114 // Face not yet cut so insert.
2115 newFaceSplitCut.insert(meshFacei, cutEdge);
2116 }
2117 else
2118 {
2119 // Face already cut. Ok if same edge.
2120 if (iter.val() != cutEdge)
2121 {
2122 return false;
2123 }
2124 }
2125 }
2126 }
2127
2128 // Is there a face on which all cuts are?
2129 label faceContainingLoop = loopFace(celli, loop);
2130
2131 if (faceContainingLoop != -1)
2132 {
2133 if (verbose_ || debug)
2134 {
2136 << "Found loop on cell " << celli << " with all points"
2137 << " on face " << faceContainingLoop << endl;
2138 }
2139
2140 //writeOBJ(".", celli, loopPoints(loop, loopWeights), labelList(0));
2141
2142 return false;
2143 }
2144
2145 // Calculate anchor points
2146 // Final success is determined by whether anchor points can be determined.
2147 return calcAnchors
2148 (
2149 celli,
2150 loop,
2151 loopPoints(loop, loopWeights),
2152 anchorPoints
2153 );
2154}
2155
2156
2157void Foam::cellCuts::setFromCellLoops()
2158{
2159 // 'Uncut' edges/vertices that are not used in loops.
2160 pointIsCut_ = false;
2161
2162 edgeIsCut_ = false;
2163
2164 faceSplitCut_.clear();
2165
2166 forAll(cellLoops_, celli)
2167 {
2168 const labelList& loop = cellLoops_[celli];
2169
2170 if (loop.size())
2171 {
2172 // Storage for cross-face cuts
2173 Map<edge> faceSplitCuts(loop.size());
2174 // Storage for points on one side of cell.
2175 labelList anchorPoints;
2176
2177 if
2178 (
2179 !validLoop
2180 (
2181 celli,
2182 loop,
2183 loopWeights(loop),
2184 faceSplitCuts,
2185 anchorPoints
2186 )
2187 )
2188 {
2189 //writeOBJ(".", celli, loopPoints(celli), anchorPoints);
2190 if (verbose_ || debug)
2191 {
2193 << "Illegal loop " << loop
2194 << " when recreating cut-addressing"
2195 << " from existing cellLoops for cell " << celli
2196 << endl;
2197 }
2198
2199 cellLoops_[celli].clear();
2200 cellAnchorPoints_[celli].clear();
2201 }
2202 else
2203 {
2204 // Copy anchor points.
2205 cellAnchorPoints_[celli].transfer(anchorPoints);
2206
2207 // Copy faceSplitCuts into overall faceSplit info.
2208 forAllConstIters(faceSplitCuts, iter)
2209 {
2210 faceSplitCut_.insert(iter.key(), iter.val());
2211 }
2212
2213 // Update edgeIsCut, pointIsCut information
2214 forAll(loop, cutI)
2215 {
2216 label cut = loop[cutI];
2217
2218 if (isEdge(cut))
2219 {
2220 edgeIsCut_[getEdge(cut)] = true;
2221 }
2222 else
2223 {
2224 pointIsCut_[getVertex(cut)] = true;
2225 }
2226 }
2227 }
2228 }
2229 }
2230
2231 // Reset edge weights
2232 forAll(edgeIsCut_, edgeI)
2233 {
2234 if (!edgeIsCut_[edgeI])
2235 {
2236 // Weight not used. Set to illegal value to make any use fall over.
2237 edgeWeight_[edgeI] = -GREAT;
2238 }
2239 }
2240}
2241
2242
2243bool Foam::cellCuts::setFromCellLoop
2244(
2245 const label celli,
2246 const labelList& loop,
2247 const scalarField& loopWeights
2248)
2249{
2250 // Update basic cut information from single cellLoop. Returns true if loop
2251 // was valid.
2252
2253 // Dump loop for debugging.
2254 if (debug)
2255 {
2256 OFstream str("last_cell.obj");
2257
2258 str<< "# edges of cell " << celli << nl;
2259
2261 (
2262 str,
2263 mesh().cells(),
2264 mesh().faces(),
2265 mesh().points(),
2266 labelList(1, celli)
2267 );
2268
2269
2270 OFstream loopStr("last_loop.obj");
2271
2272 loopStr<< "# looppoints for cell " << celli << nl;
2273
2274 pointField pointsOfLoop = loopPoints(loop, loopWeights);
2275
2276 forAll(pointsOfLoop, i)
2277 {
2278 meshTools::writeOBJ(loopStr, pointsOfLoop[i]);
2279 }
2280
2281 str << 'l';
2282
2283 forAll(pointsOfLoop, i)
2284 {
2285 str << ' ' << i + 1;
2286 }
2287 str << ' ' << 1 << nl;
2288 }
2289
2290 bool okLoop = false;
2291
2292 if (validEdgeLoop(loop, loopWeights))
2293 {
2294 // Storage for cuts across face
2295 Map<edge> faceSplitCuts(loop.size());
2296 // Storage for points on one side of cell
2297 labelList anchorPoints;
2298
2299 okLoop =
2300 validLoop(celli, loop, loopWeights, faceSplitCuts, anchorPoints);
2301
2302 if (okLoop)
2303 {
2304 // Valid loop. Copy cellLoops and anchorPoints
2305 cellLoops_[celli] = loop;
2306 cellAnchorPoints_[celli].transfer(anchorPoints);
2307
2308 // Copy split cuts
2309 forAllConstIters(faceSplitCuts, iter)
2310 {
2311 faceSplitCut_.insert(iter.key(), iter.val());
2312 }
2313
2314
2315 // Update edgeIsCut, pointIsCut information
2316 forAll(loop, cutI)
2317 {
2318 const label cut = loop[cutI];
2319
2320 if (isEdge(cut))
2321 {
2322 const label edgeI = getEdge(cut);
2323
2324 edgeIsCut_[edgeI] = true;
2325
2326 edgeWeight_[edgeI] = loopWeights[cutI];
2327 }
2328 else
2329 {
2330 const label vertI = getVertex(cut);
2331
2332 pointIsCut_[vertI] = true;
2333 }
2334 }
2335 }
2336 }
2337
2338 return okLoop;
2339}
2340
2341
2342void Foam::cellCuts::setFromCellLoops
2343(
2344 const labelList& cellLabels,
2345 const labelListList& cellLoops,
2346 const List<scalarField>& cellLoopWeights
2347)
2348{
2349 // 'Uncut' edges/vertices that are not used in loops.
2350 pointIsCut_ = false;
2351
2352 edgeIsCut_ = false;
2353
2354 forAll(cellLabels, cellLabelI)
2355 {
2356 const label celli = cellLabels[cellLabelI];
2357
2358 const labelList& loop = cellLoops[cellLabelI];
2359
2360 if (loop.size())
2361 {
2362 const scalarField& loopWeights = cellLoopWeights[cellLabelI];
2363
2364 if (setFromCellLoop(celli, loop, loopWeights))
2365 {
2366 // Valid loop. Call above will have updated all already.
2367 }
2368 else
2369 {
2370 // Clear cellLoops
2371 cellLoops_[celli].clear();
2372 }
2373 }
2374 }
2375}
2376
2377
2378void Foam::cellCuts::setFromCellCutter
2379(
2380 const cellLooper& cellCutter,
2381 const List<refineCell>& refCells
2382)
2383{
2384 // 'Uncut' edges/vertices that are not used in loops.
2385 pointIsCut_ = false;
2386
2387 edgeIsCut_ = false;
2388
2389 // storage for loop of cuts (cut vertices and/or cut edges)
2390 labelList cellLoop;
2391 scalarField cellLoopWeights;
2392
2393 // For debugging purposes
2394 DynamicList<label> invalidCutCells(2);
2395 DynamicList<labelList> invalidCutLoops(2);
2396 DynamicList<scalarField> invalidCutLoopWeights(2);
2397
2398
2399 forAll(refCells, refCelli)
2400 {
2401 const refineCell& refCell = refCells[refCelli];
2402
2403 const label celli = refCell.cellNo();
2404
2405 const vector& refDir = refCell.direction();
2406
2407
2408 // Cut cell. Determines cellLoop and cellLoopWeights
2409 bool goodCut =
2410 cellCutter.cut
2411 (
2412 refDir,
2413 celli,
2414
2415 pointIsCut_,
2416 edgeIsCut_,
2417 edgeWeight_,
2418
2419 cellLoop,
2420 cellLoopWeights
2421 );
2422
2423 // Check whether edge refinement is on a per face basis compatible with
2424 // current pattern.
2425 if (goodCut)
2426 {
2427 if (setFromCellLoop(celli, cellLoop, cellLoopWeights))
2428 {
2429 // Valid loop. Will have updated all info already.
2430 }
2431 else
2432 {
2433 cellLoops_[celli].clear();
2434
2436 << "Found loop on cell " << celli
2437 << " that resulted in an unexpected bad cut." << nl
2438 << " Suggestions:" << nl
2439 << " - Turn on the debug switch for 'cellCuts' to get"
2440 << " geometry files that identify this cell." << nl
2441 << " - Also keep in mind to check the defined"
2442 << " reference directions, as these are most likely the"
2443 << " origin of the problem."
2444 << nl << endl;
2445
2446 // Discarded by validLoop
2447 if (debug)
2448 {
2449 invalidCutCells.append(celli);
2450 invalidCutLoops.append(cellLoop);
2451 invalidCutLoopWeights.append(cellLoopWeights);
2452 }
2453 }
2454 }
2455 else
2456 {
2457 // Clear cellLoops
2458 cellLoops_[celli].clear();
2459 }
2460 }
2461
2462 if (debug && invalidCutCells.size())
2463 {
2464 invalidCutCells.shrink();
2465 invalidCutLoops.shrink();
2466 invalidCutLoopWeights.shrink();
2467
2468 fileName cutsFile("invalidLoopCells.obj");
2469
2470 Pout<< "cellCuts : writing inValidLoops cells to " << cutsFile << endl;
2471
2472 OFstream cutsStream(cutsFile);
2473
2475 (
2476 cutsStream,
2477 mesh().cells(),
2478 mesh().faces(),
2479 mesh().points(),
2480 invalidCutCells
2481 );
2482
2483 fileName loopsFile("invalidLoops.obj");
2484
2485 Pout<< "cellCuts : writing inValidLoops loops to " << loopsFile << endl;
2486
2487 OFstream loopsStream(loopsFile);
2488
2489 label vertI = 0;
2490
2491 forAll(invalidCutLoops, i)
2492 {
2493 writeOBJ
2494 (
2495 loopsStream,
2496 loopPoints(invalidCutLoops[i], invalidCutLoopWeights[i]),
2497 vertI
2498 );
2499 }
2500 }
2501}
2502
2503
2504void Foam::cellCuts::setFromCellCutter
2505(
2506 const cellLooper& cellCutter,
2507 const labelList& cellLabels,
2508 const List<plane>& cellCutPlanes
2509)
2510{
2511 // 'Uncut' edges/vertices that are not used in loops.
2512 pointIsCut_ = false;
2513
2514 edgeIsCut_ = false;
2515
2516 // storage for loop of cuts (cut vertices and/or cut edges)
2517 labelList cellLoop;
2518 scalarField cellLoopWeights;
2519
2520 // For debugging purposes
2521 DynamicList<label> invalidCutCells(2);
2522 DynamicList<labelList> invalidCutLoops(2);
2523 DynamicList<scalarField> invalidCutLoopWeights(2);
2524
2525
2526 forAll(cellLabels, i)
2527 {
2528 const label celli = cellLabels[i];
2529
2530 // Cut cell. Determines cellLoop and cellLoopWeights
2531 bool goodCut =
2532 cellCutter.cut
2533 (
2534 cellCutPlanes[i],
2535 celli,
2536
2537 pointIsCut_,
2538 edgeIsCut_,
2539 edgeWeight_,
2540
2541 cellLoop,
2542 cellLoopWeights
2543 );
2544
2545 // Check whether edge refinement is on a per face basis compatible with
2546 // current pattern.
2547 if (goodCut)
2548 {
2549 if (setFromCellLoop(celli, cellLoop, cellLoopWeights))
2550 {
2551 // Valid loop. Will have updated all info already.
2552 }
2553 else
2554 {
2555 cellLoops_[celli].clear();
2556
2557 // Discarded by validLoop
2558 if (debug)
2559 {
2560 invalidCutCells.append(celli);
2561 invalidCutLoops.append(cellLoop);
2562 invalidCutLoopWeights.append(cellLoopWeights);
2563 }
2564 }
2565 }
2566 else
2567 {
2568 // Clear cellLoops
2569 cellLoops_[celli].clear();
2570 }
2571 }
2572
2573 if (debug && invalidCutCells.size())
2574 {
2575 invalidCutCells.shrink();
2576 invalidCutLoops.shrink();
2577 invalidCutLoopWeights.shrink();
2578
2579 fileName cutsFile("invalidLoopCells.obj");
2580
2581 Pout<< "cellCuts : writing inValidLoops cells to " << cutsFile << endl;
2582
2583 OFstream cutsStream(cutsFile);
2584
2586 (
2587 cutsStream,
2588 mesh().cells(),
2589 mesh().faces(),
2590 mesh().points(),
2591 invalidCutCells
2592 );
2593
2594 fileName loopsFile("invalidLoops.obj");
2595
2596 Pout<< "cellCuts : writing inValidLoops loops to " << loopsFile << endl;
2597
2598 OFstream loopsStream(loopsFile);
2599
2600 label vertI = 0;
2601
2602 forAll(invalidCutLoops, i)
2603 {
2604 writeOBJ
2605 (
2606 loopsStream,
2607 loopPoints(invalidCutLoops[i], invalidCutLoopWeights[i]),
2608 vertI
2609 );
2610 }
2611 }
2612}
2613
2614
2615void Foam::cellCuts::orientPlanesAndLoops()
2616{
2617 // Determine anchorPoints if not yet done by validLoop.
2618 forAll(cellLoops_, celli)
2619 {
2620 const labelList& loop = cellLoops_[celli];
2621
2622 if (loop.size() && cellAnchorPoints_[celli].empty())
2623 {
2624 // Leave anchor points empty if illegal loop.
2625 calcAnchors
2626 (
2627 celli,
2628 loop,
2629 loopPoints(celli),
2630 cellAnchorPoints_[celli]
2631 );
2632 }
2633 }
2634
2635 if (debug & 2)
2636 {
2637 Pout<< "cellAnchorPoints:" << endl;
2638 }
2639 forAll(cellAnchorPoints_, celli)
2640 {
2641 if (cellLoops_[celli].size())
2642 {
2643 if (cellAnchorPoints_[celli].empty())
2644 {
2646 << "No anchor points for cut cell " << celli << endl
2647 << "cellLoop:" << cellLoops_[celli] << abort(FatalError);
2648 }
2649
2650 if (debug & 2)
2651 {
2652 Pout<< " cell:" << celli << " anchored at "
2653 << cellAnchorPoints_[celli] << endl;
2654 }
2655 }
2656 }
2657
2658 // Calculate number of valid cellLoops
2659 nLoops_ = 0;
2660
2661 forAll(cellLoops_, celli)
2662 {
2663 if (cellLoops_[celli].size())
2664 {
2665 nLoops_++;
2666 }
2667 }
2668}
2669
2670
2671void Foam::cellCuts::calcLoopsAndAddressing(const labelList& cutCells)
2672{
2673 // Sanity check on weights
2674 forAll(edgeIsCut_, edgeI)
2675 {
2676 if (edgeIsCut_[edgeI])
2677 {
2678 scalar weight = edgeWeight_[edgeI];
2679
2680 if (weight < 0 || weight > 1)
2681 {
2683 << "Weight out of range [0,1]. Edge " << edgeI
2684 << " verts:" << mesh().edges()[edgeI]
2685 << " weight:" << weight << abort(FatalError);
2686 }
2687 }
2688 else
2689 {
2690 // Weight not used. Set to illegal value to make any use fall over.
2691 edgeWeight_[edgeI] = -GREAT;
2692 }
2693 }
2694
2695
2696 // Calculate faces that split cells in two
2697 calcCellLoops(cutCells);
2698
2699 if (debug & 2)
2700 {
2701 Pout<< "-- cellLoops --" << endl;
2702 forAll(cellLoops_, celli)
2703 {
2704 const labelList& loop = cellLoops_[celli];
2705
2706 if (loop.size())
2707 {
2708 Pout<< "cell:" << celli << " ";
2709 writeCuts(Pout, loop, loopWeights(loop));
2710 Pout<< endl;
2711 }
2712 }
2713 }
2714
2715 // Redo basic cut information (pointIsCut, edgeIsCut, faceSplitCut)
2716 // using cellLoop only.
2717 setFromCellLoops();
2718}
2719
2720
2721void Foam::cellCuts::check() const
2722{
2723 // Check weights for unsnapped values
2724 forAll(edgeIsCut_, edgeI)
2725 {
2726 if (edgeIsCut_[edgeI])
2727 {
2728 if
2729 (
2730 edgeWeight_[edgeI] < geomCellLooper::snapTol()
2731 || edgeWeight_[edgeI] > (1 - geomCellLooper::snapTol())
2732 )
2733 {
2734 // Should have been snapped.
2736 << "edge:" << edgeI << " vertices:"
2737 << mesh().edges()[edgeI]
2738 << " weight:" << edgeWeight_[edgeI] << " should have been"
2739 << " snapped to one of its endpoints"
2740 << endl;
2741 }
2742 }
2743 else
2744 {
2745 if (edgeWeight_[edgeI] > - 1)
2746 {
2748 << "edge:" << edgeI << " vertices:"
2749 << mesh().edges()[edgeI]
2750 << " weight:" << edgeWeight_[edgeI] << " is not cut but"
2751 << " its weight is not marked invalid"
2752 << abort(FatalError);
2753 }
2754 }
2755 }
2756
2757 // Check that all elements of cellloop are registered
2758 forAll(cellLoops_, celli)
2759 {
2760 const labelList& loop = cellLoops_[celli];
2761
2762 forAll(loop, i)
2763 {
2764 label cut = loop[i];
2765
2766 if
2767 (
2768 (isEdge(cut) && !edgeIsCut_[getEdge(cut)])
2769 || (!isEdge(cut) && !pointIsCut_[getVertex(cut)])
2770 )
2771 {
2772 labelList cuts(1, cut);
2773 writeCuts(Pout, cuts, loopWeights(cuts));
2774
2776 << "cell:" << celli << " loop:"
2777 << loop
2778 << " cut:" << cut << " is not marked as cut"
2779 << abort(FatalError);
2780 }
2781 }
2782 }
2783
2784 // Check that no elements of cell loop are anchor point.
2785 forAll(cellLoops_, celli)
2786 {
2787 const labelList& anchors = cellAnchorPoints_[celli];
2788
2789 const labelList& loop = cellLoops_[celli];
2790
2791 if (loop.size() && anchors.empty())
2792 {
2794 << "cell:" << celli << " loop:" << loop
2795 << " has no anchor points"
2796 << abort(FatalError);
2797 }
2798
2799
2800 forAll(loop, i)
2801 {
2802 label cut = loop[i];
2803
2804 if
2805 (
2806 !isEdge(cut)
2807 && anchors.found(getVertex(cut))
2808 )
2809 {
2811 << "cell:" << celli << " loop:" << loop
2812 << " anchor points:" << anchors
2813 << " anchor:" << getVertex(cut) << " is part of loop"
2814 << abort(FatalError);
2815 }
2816 }
2817 }
2818
2819
2820 // Check that cut faces have a neighbour that is cut.
2821 boolList nbrCellIsCut;
2822 {
2823 boolList cellIsCut(mesh().nCells(), false);
2824 forAll(cellLoops_, celli)
2825 {
2826 cellIsCut[celli] = cellLoops_[celli].size();
2827 }
2828 syncTools::swapBoundaryCellList(mesh(), cellIsCut, nbrCellIsCut);
2829 }
2830
2831 forAllConstIters(faceSplitCut_, iter)
2832 {
2833 const label facei = iter.key();
2834
2835 if (mesh().isInternalFace(facei))
2836 {
2837 label own = mesh().faceOwner()[facei];
2838 label nei = mesh().faceNeighbour()[facei];
2839
2840 if (cellLoops_[own].empty() && cellLoops_[nei].empty())
2841 {
2843 << "Internal face:" << facei << " cut by " << iter.val()
2844 << " has owner:" << own
2845 << " and neighbour:" << nei
2846 << " that are both uncut"
2847 << abort(FatalError);
2848 }
2849 }
2850 else
2851 {
2852 label bFacei = facei - mesh().nInternalFaces();
2853
2854 label own = mesh().faceOwner()[facei];
2855
2856 if (cellLoops_[own].empty() && !nbrCellIsCut[bFacei])
2857 {
2859 << "Boundary face:" << facei << " cut by " << iter.val()
2860 << " has owner:" << own
2861 << " that is uncut"
2862 << abort(FatalError);
2863 }
2864 }
2865 }
2866}
2867
2868
2869// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
2870
2872(
2873 const polyMesh& mesh,
2874 const labelList& cutCells,
2875 const labelList& meshVerts,
2876 const labelList& meshEdges,
2877 const scalarField& meshEdgeWeights,
2878 const bool verbose
2879)
2880:
2882 verbose_(verbose),
2883 pointIsCut_(expand(mesh.nPoints(), meshVerts)),
2884 edgeIsCut_(expand(mesh.nEdges(), meshEdges)),
2885 edgeWeight_(expand(mesh.nEdges(), meshEdges, meshEdgeWeights)),
2886 faceSplitCut_(cutCells.size()),
2887 cellLoops_(mesh.nCells()),
2888 nLoops_(-1),
2889 cellAnchorPoints_(mesh.nCells())
2890{
2891 if (debug)
2892 {
2893 Pout<< "cellCuts : constructor from cut verts and edges" << endl;
2894 }
2895
2896 calcLoopsAndAddressing(cutCells);
2897
2898 // Calculate planes and flip cellLoops if necessary
2899 orientPlanesAndLoops();
2900
2901 if (debug)
2902 {
2903 check();
2904 }
2905
2906 clearOut();
2907
2908 if (debug)
2909 {
2910 Pout<< "cellCuts : leaving constructor from cut verts and edges"
2911 << endl;
2912 }
2913}
2914
2915
2917(
2918 const polyMesh& mesh,
2919 const labelList& meshVerts,
2920 const labelList& meshEdges,
2921 const scalarField& meshEdgeWeights,
2922 const bool verbose
2923)
2924:
2926 verbose_(verbose),
2927 pointIsCut_(expand(mesh.nPoints(), meshVerts)),
2928 edgeIsCut_(expand(mesh.nEdges(), meshEdges)),
2929 edgeWeight_(expand(mesh.nEdges(), meshEdges, meshEdgeWeights)),
2930 faceSplitCut_(mesh.nFaces()/10 + 1),
2931 cellLoops_(mesh.nCells()),
2932 nLoops_(-1),
2933 cellAnchorPoints_(mesh.nCells())
2934{
2935 // Construct from pattern of cuts. Finds out itself which cells are cut.
2936 // (can go wrong if e.g. all neighbours of cell are refined)
2937
2938 if (debug)
2939 {
2940 Pout<< "cellCuts : constructor from cellLoops" << endl;
2941 }
2942
2943 calcLoopsAndAddressing(identity(mesh.nCells()));
2944
2945 // Adds cuts on other side of coupled boundaries
2946 syncProc();
2947
2948 // Calculate planes and flip cellLoops if necessary
2949 orientPlanesAndLoops();
2950
2951 if (debug)
2952 {
2953 check();
2954 }
2955
2956 clearOut();
2957
2958 if (debug)
2959 {
2960 Pout<< "cellCuts : leaving constructor from cellLoops" << endl;
2961 }
2962}
2963
2964
2966(
2967 const polyMesh& mesh,
2968 const labelList& cellLabels,
2969 const labelListList& cellLoops,
2970 const List<scalarField>& cellEdgeWeights,
2971 const bool verbose
2972)
2973:
2975 verbose_(verbose),
2976 pointIsCut_(mesh.nPoints(), false),
2977 edgeIsCut_(mesh.nEdges(), false),
2978 edgeWeight_(mesh.nEdges(), -GREAT),
2979 faceSplitCut_(cellLabels.size()),
2980 cellLoops_(mesh.nCells()),
2981 nLoops_(-1),
2982 cellAnchorPoints_(mesh.nCells())
2983{
2984 if (debug)
2985 {
2986 Pout<< "cellCuts : constructor from cellLoops" << endl;
2987 }
2988
2989 // Update pointIsCut, edgeIsCut, faceSplitCut from cell loops.
2990 // Makes sure cuts are consistent
2991 setFromCellLoops(cellLabels, cellLoops, cellEdgeWeights);
2992
2993 // Adds cuts on other side of coupled boundaries
2994 syncProc();
2995
2996 // Calculate planes and flip cellLoops if necessary
2997 orientPlanesAndLoops();
2998
2999 if (debug)
3000 {
3001 check();
3002 }
3003
3004 clearOut();
3005
3006 if (debug)
3007 {
3008 Pout<< "cellCuts : leaving constructor from cellLoops" << endl;
3009 }
3010}
3011
3012
3014(
3015 const polyMesh& mesh,
3016 const cellLooper& cellCutter,
3017 const List<refineCell>& refCells,
3018 const bool verbose
3019)
3020:
3022 verbose_(verbose),
3023 pointIsCut_(mesh.nPoints(), false),
3024 edgeIsCut_(mesh.nEdges(), false),
3025 edgeWeight_(mesh.nEdges(), -GREAT),
3026 faceSplitCut_(refCells.size()),
3027 cellLoops_(mesh.nCells()),
3028 nLoops_(-1),
3029 cellAnchorPoints_(mesh.nCells())
3030{
3031 if (debug)
3032 {
3033 Pout<< "cellCuts : constructor from cellCutter" << endl;
3034 }
3035
3036 // Update pointIsCut, edgeIsCut, faceSplitCut from cell loops.
3037 // Makes sure cuts are consistent
3038 setFromCellCutter(cellCutter, refCells);
3039
3040 // Adds cuts on other side of coupled boundaries
3041 syncProc();
3042
3043 // Calculate planes and flip cellLoops if necessary
3044 orientPlanesAndLoops();
3045
3046 if (debug)
3047 {
3048 check();
3049 }
3050
3051 clearOut();
3052
3053 if (debug)
3054 {
3055 Pout<< "cellCuts : leaving constructor from cellCutter" << endl;
3056 }
3057}
3058
3059
3061(
3062 const polyMesh& mesh,
3063 const cellLooper& cellCutter,
3064 const labelList& cellLabels,
3065 const List<plane>& cutPlanes,
3066 const bool verbose
3067)
3068:
3070 verbose_(verbose),
3071 pointIsCut_(mesh.nPoints(), false),
3072 edgeIsCut_(mesh.nEdges(), false),
3073 edgeWeight_(mesh.nEdges(), -GREAT),
3074 faceSplitCut_(cellLabels.size()),
3075 cellLoops_(mesh.nCells()),
3076 nLoops_(-1),
3077 cellAnchorPoints_(mesh.nCells())
3078{
3079 if (debug)
3080 {
3081 Pout<< "cellCuts : constructor from cellCutter with prescribed plane"
3082 << endl;
3083 }
3084
3085 // Update pointIsCut, edgeIsCut, faceSplitCut from cell loops.
3086 // Makes sure cuts are consistent
3087 setFromCellCutter(cellCutter, cellLabels, cutPlanes);
3088
3089 // Adds cuts on other side of coupled boundaries
3090 syncProc();
3091
3092 // Calculate planes and flip cellLoops if necessary
3093 orientPlanesAndLoops();
3094
3095 if (debug)
3096 {
3097 check();
3098 }
3099
3100 clearOut();
3101
3102 if (debug)
3103 {
3104 Pout<< "cellCuts : leaving constructor from cellCutter with prescribed"
3105 << " plane" << endl;
3106 }
3107}
3108
3109
3111(
3112 const polyMesh& mesh,
3113 const boolList& pointIsCut,
3114 const boolList& edgeIsCut,
3115 const scalarField& edgeWeight,
3116 const Map<edge>& faceSplitCut,
3117 const labelListList& cellLoops,
3118 const label nLoops,
3119 const labelListList& cellAnchorPoints,
3120 const bool verbose
3121)
3122:
3124 verbose_(verbose),
3125 pointIsCut_(pointIsCut),
3126 edgeIsCut_(edgeIsCut),
3127 edgeWeight_(edgeWeight),
3128 faceSplitCut_(faceSplitCut),
3129 cellLoops_(cellLoops),
3130 nLoops_(nLoops),
3131 cellAnchorPoints_(cellAnchorPoints)
3132{
3133 if (debug)
3134 {
3135 Pout<< "cellCuts : constructor from components" << endl;
3136 Pout<< "cellCuts : leaving constructor from components" << endl;
3137 }
3138}
3139
3140
3141// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
3142
3144{
3145 faceCutsPtr_.reset(nullptr);
3146}
3147
3148
3149Foam::pointField Foam::cellCuts::loopPoints(const label celli) const
3150{
3151 const labelList& loop = cellLoops_[celli];
3152
3153 pointField loopPts(loop.size());
3154
3155 forAll(loop, fp)
3156 {
3157 const label cut = loop[fp];
3158
3159 if (isEdge(cut))
3160 {
3161 label edgeI = getEdge(cut);
3162
3163 loopPts[fp] = coord(cut, edgeWeight_[edgeI]);
3164 }
3165 else
3166 {
3167 loopPts[fp] = coord(cut, -GREAT);
3168 }
3169 }
3170 return loopPts;
3171}
3172
3173
3174void Foam::cellCuts::flip(const label celli)
3175{
3176 labelList& loop = cellLoops_[celli];
3177
3178 reverse(loop);
3179
3180 // Reverse anchor point set.
3181 cellAnchorPoints_[celli] =
3182 nonAnchorPoints
3183 (
3184 mesh().cellPoints()[celli],
3185 cellAnchorPoints_[celli],
3186 loop
3187 );
3188}
3189
3190
3191void Foam::cellCuts::flipLoopOnly(const label celli)
3192{
3193 labelList& loop = cellLoops_[celli];
3194
3195 reverse(loop);
3196}
3197
3198
3199void Foam::cellCuts::writeOBJ
3200(
3201 Ostream& os,
3202 const pointField& loopPts,
3203 label& vertI
3204) const
3205{
3206 label startVertI = vertI;
3207
3208 forAll(loopPts, fp)
3209 {
3210 const point& pt = loopPts[fp];
3211
3212 os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
3213
3214 vertI++;
3215 }
3216
3217 os << 'f';
3218 forAll(loopPts, fp)
3219 {
3220 os << ' ' << startVertI + fp + 1;
3221 }
3222 os << endl;
3223}
3224
3225
3226void Foam::cellCuts::writeOBJ(Ostream& os) const
3227{
3228 label vertI = 0;
3229
3230 forAll(cellLoops_, celli)
3231 {
3232 if (cellLoops_[celli].size())
3233 {
3234 writeOBJ(os, loopPoints(celli), vertI);
3235 }
3236 }
3237}
3238
3239
3240void Foam::cellCuts::writeCellOBJ(const fileName& dir, const label celli) const
3241{
3242 const labelList& anchors = cellAnchorPoints_[celli];
3243
3244 writeOBJ(dir, celli, loopPoints(celli), anchors);
3245}
3246
3247
3248// ************************************************************************* //
Various functions to operate on Lists.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
bool found(const T &val, label pos=0) const
True if the value if found in the list.
Definition: UListI.H:265
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
label rcIndex(const label i) const noexcept
Definition: UListI.H:67
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label fcIndex(const label i) const noexcept
Definition: UListI.H:60
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
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
void writeCellOBJ(const fileName &dir, const label celli) const
debugging:Write edges of cell and loop
Definition: cellCuts.C:3240
void flipLoopOnly(const label celli)
Flip loop for celli. Does not update anchors. Use with care.
Definition: cellCuts.C:3191
labelList nonAnchorPoints(const labelList &cellPoints, const labelList &anchorPoints, const labelList &loop) const
Invert anchor point selection.
Definition: cellCuts.C:1330
void clearOut()
Clear out demand-driven storage.
Definition: cellCuts.C:3143
Abstract base class. Concrete implementations know how to cut a cell (i.e. determine a loop around th...
Definition: cellLooper.H:75
Combines edge or vertex in single label. Used to specify cuts across cell circumference.
Definition: edgeVertex.H:56
const polyMesh & mesh() const
Definition: edgeVertex.H:101
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A class for handling file names.
Definition: fileName.H:76
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
static scalar snapTol()
edge cmptType
Component type.
Definition: cellCuts.C:56
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
void flip()
Reverse the orientation.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
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
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 labelListList & cellPoints() const
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
const labelListList & pointFaces() const
const labelListList & edgeFaces() const
const labelListList & faceEdges() const
const cellList & cells() const
Patchify triangles based on orientation w.r.t other (triangulated or triangulatable) surfaces.
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
static void syncEdgeList(const polyMesh &mesh, List< T > &edgeValues, const CombineOp &cop, const T &nullValue, const TransformOp &top, const FlipOp &fop)
Synchronize values on all mesh edges.
static void syncBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const TransformOp &top, const bool parRun=Pstream::parRun())
Synchronize values on boundary faces only.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
Dummy transform to be used with syncTools.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
const labelList nEdges(UPstream::listGatherValues< label >(aMesh.nEdges()))
const pointField & points
label nPoints
const cellShapeList & cells
#define WarningInFunction
Report a warning using Foam::Warning.
bool faceOnCell(const primitiveMesh &mesh, const label celli, const label facei)
Is face used by cell.
Definition: meshTools.C:330
bool edgeOnCell(const primitiveMesh &mesh, const label celli, const label edgeI)
Is edge used by cell.
Definition: meshTools.C:308
label walkFace(const primitiveMesh &mesh, const label facei, const label startEdgeI, const label startVertI, const label nEdges)
Returns label of edge nEdges away from startEdge (in the direction.
Definition: meshTools.C:603
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
label otherFace(const primitiveMesh &mesh, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:555
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:479
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Namespace for OpenFOAM.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
List< label > labelList
A List of labels.
Definition: List.H:66
constexpr label labelMin
Definition: label.H:60
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition: UListI.H:449
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:144
List< bool > boolList
A List of bools.
Definition: List.H:64
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
error FatalError
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
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)
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