intersectedSurface.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) 2015-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 "intersectedSurface.H"
30#include "surfaceIntersection.H"
31#include "faceList.H"
32#include "faceTriangulation.H"
33#include "treeBoundBox.H"
34#include "OFstream.H"
35#include "error.H"
36#include "meshTools.H"
37#include "edgeSurface.H"
38#include "DynamicList.H"
39#include "transform.H"
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43namespace Foam
44{
46}
47
48
49const Foam::label Foam::intersectedSurface::UNVISITED = 0;
50const Foam::label Foam::intersectedSurface::STARTTOEND = 1;
51const Foam::label Foam::intersectedSurface::ENDTOSTART = 2;
52const Foam::label Foam::intersectedSurface::BOTH = STARTTOEND | ENDTOSTART;
53
54
55// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56
57void Foam::intersectedSurface::writeOBJ
58(
59 const pointField& points,
60 const edgeList& edges,
61 Ostream& os
62)
63{
64 for (const point& pt : points)
65 {
66 os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
67 }
68
69 for (const edge& e : edges)
70 {
71 os << "l " << e.start()+1 << ' ' << e.end()+1 << nl;
72 }
73}
74
75
76void Foam::intersectedSurface::writeOBJ
77(
78 const pointField& points,
79 const edgeList& edges,
80 const labelList& faceEdges,
81 Ostream& os
82)
83{
84 for (const point& pt : points)
85 {
86 os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
87 }
88 for (const label edgei : faceEdges)
89 {
90 const edge& e = edges[edgei];
91
92 os << "l " << e.start()+1 << ' ' << e.end()+1 << nl;
93 }
94}
95
96
97void Foam::intersectedSurface::writeLocalOBJ
98(
99 const pointField& points,
100 const edgeList& edges,
101 const labelList& faceEdges,
102 const fileName& fName
103)
104{
105 OFstream os(fName);
106
107 labelList pointMap(points.size(), -1);
108
109 label maxVerti = 0;
110
111 for (const label edgei : faceEdges)
112 {
113 const edge& e = edges[edgei];
114
115 forAll(e, i)
116 {
117 const label pointi = e[i];
118
119 if (pointMap[pointi] == -1)
120 {
121 const point& pt = points[pointi];
122
123 os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
124
125 pointMap[pointi] = maxVerti++;
126 }
127 }
128 }
129
130 for (const label edgei : faceEdges)
131 {
132 const edge& e = edges[edgei];
133
134 os << "l " << pointMap[e.start()]+1 << ' ' << pointMap[e.end()]+1
135 << nl;
136 }
137}
138
139
140void Foam::intersectedSurface::writeOBJ
141(
142 const pointField& points,
143 const face& f,
144 Ostream& os
145)
146{
147 for (const point& pt : points)
148 {
149 os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
150 }
151
152 os << 'f';
153
154 for (const label pointi : f)
155 {
156 os << ' ' << pointi+1;
157 }
158 os << nl;
159}
160
161
162void Foam::intersectedSurface::printVisit
163(
164 const edgeList& edges,
165 const labelList& edgeLabels,
166 const Map<label>& visited
167)
168{
169 Pout<< "Visited:" << nl;
170 for (const label edgei : edgeLabels)
171 {
172 const edge& e = edges[edgei];
173
174 label stat = visited[edgei];
175
176 if (stat == UNVISITED)
177 {
178 Pout<< " edge:" << edgei << " verts:" << e
179 << " unvisited" << nl;
180 }
181 else if (stat == STARTTOEND)
182 {
183 Pout<< " edge:" << edgei << " from " << e[0]
184 << " to " << e[1] << nl;
185 }
186 else if (stat == ENDTOSTART)
187 {
188 Pout<< " edge:" << edgei << " from " << e[1]
189 << " to " << e[0] << nl;
190 }
191 else
192 {
193 Pout<< " edge:" << edgei << " both " << e
194 << nl;
195 }
196 }
197 Pout<< endl;
198}
199
200
201bool Foam::intersectedSurface::sameEdgeOrder
202(
203 const labelledTri& fA,
204 const labelledTri& fB
205)
206{
207 forAll(fA, fpA)
208 {
209 label fpB = fB.find(fA[fpA]);
210
211 if (fpB != -1)
212 {
213 // Get prev/next vertex on fA
214 label vA1 = fA[fA.fcIndex(fpA)];
215 label vAMin1 = fA[fA.rcIndex(fpA)];
216
217 // Get prev/next vertex on fB
218 label vB1 = fB[fB.fcIndex(fpB)];
219 label vBMin1 = fB[fB.rcIndex(fpB)];
220
221 if (vA1 == vB1 || vAMin1 == vBMin1)
222 {
223 return true;
224 }
225 else if (vA1 == vBMin1 || vAMin1 == vB1)
226 {
227 // shared vertices in opposite order.
228 return false;
229 }
230 else
231 {
233 << "Triangle:" << fA << " and triangle:" << fB
234 << " share a point but not an edge"
235 << abort(FatalError);
236 }
237 }
238 }
239
241 << "Triangle:" << fA << " and triangle:" << fB
242 << " do not share an edge"
243 << abort(FatalError);
244
245 return false;
246}
247
248
249void Foam::intersectedSurface::incCount
250(
251 Map<label>& visited,
252 const label key,
253 const label offset
254)
255{
256 visited(key, 0) += offset;
257}
258
259
261Foam::intersectedSurface::calcPointEdgeAddressing
262(
263 const edgeSurface& eSurf,
264 const label facei
265)
266{
267 const pointField& points = eSurf.points();
268 const edgeList& edges = eSurf.edges();
269
270 const labelList& fEdges = eSurf.faceEdges()[facei];
271
272 Map<DynamicList<label>> facePointEdges(4*fEdges.size());
273
274 for (const label edgei : fEdges)
275 {
276 const edge& e = edges[edgei];
277
278 // Add e.start to point-edges
279 facePointEdges(e.start()).append(edgei);
280
281 // Add e.end to point-edges
282 facePointEdges(e.end()).append(edgei);
283 }
284
285 // Shrink it
286 forAllIters(facePointEdges, iter)
287 {
288 iter.val().shrink();
289
290 // Check on dangling points.
291 if (iter.val().empty())
292 {
294 << "Point:" << iter.key() << " used by too few edges:"
295 << iter.val() << abort(FatalError);
296 }
297 }
298
299 if (debug & 2)
300 {
301 // Print facePointEdges
302 Pout<< "calcPointEdgeAddressing: face consisting of edges:" << nl;
303 for (const label edgei : fEdges)
304 {
305 const edge& e = edges[edgei];
306 Pout<< " " << edgei << ' ' << e
307 << points[e.start()]
308 << points[e.end()] << nl;
309 }
310
311 Pout<< " Constructed point-edge addressing:" << nl;
312 forAllConstIters(facePointEdges, iter)
313 {
314 Pout<< " vertex " << iter.key() << " is connected to edges "
315 << iter.val() << nl;
316 }
317 Pout<< endl;
318 }
319
320 return facePointEdges;
321}
322
323
324Foam::label Foam::intersectedSurface::nextEdge
325(
326 const edgeSurface& eSurf,
327 const Map<label>& visited,
328 const label facei,
329 const vector& n,
330 const Map<DynamicList<label>>& facePointEdges,
331 const label prevEdgei,
332 const label prevVerti
333)
334{
335 const pointField& points = eSurf.points();
336 const edgeList& edges = eSurf.edges();
337 const labelList& fEdges = eSurf.faceEdges()[facei];
338
339
340 // Edges connected to prevVerti
341 const DynamicList<label>& connectedEdges = facePointEdges[prevVerti];
342
343 if (connectedEdges.size() <= 1)
344 {
345 // Problem. Point not connected.
346 {
347 Pout<< "Writing face:" << facei << " to face.obj" << endl;
348 OFstream str("face.obj");
349 writeOBJ(points, edges, fEdges, str);
350
351 Pout<< "Writing connectedEdges edges to faceEdges.obj" << endl;
352 writeLocalOBJ(points, edges, connectedEdges, "faceEdges.obj");
353 }
354
356 << "Problem: prevVertI:" << prevVerti << " on edge " << prevEdgei
357 << " has less than 2 connected edges."
358 << " connectedEdges:" << connectedEdges << abort(FatalError);
359
360 return -1;
361 }
362
363 if (connectedEdges.size() == 2)
364 {
365 // Simple case. Take other edge
366 if (connectedEdges[0] == prevEdgei)
367 {
368 if (debug & 2)
369 {
370 Pout<< "Stepped from edge:" << edges[prevEdgei]
371 << " to edge:" << edges[connectedEdges[1]]
372 << " chosen from candidates:" << connectedEdges << endl;
373 }
374 return connectedEdges[1];
375 }
376 else
377 {
378 if (debug & 2)
379 {
380 Pout<< "Stepped from edge:" << edges[prevEdgei]
381 << " to edge:" << edges[connectedEdges[0]]
382 << " chosen from candidates:" << connectedEdges << endl;
383 }
384 return connectedEdges[0];
385 }
386 }
387
388
389 // Multiple choices. Look at angle between edges.
390
391 const edge& prevE = edges[prevEdgei];
392
393 // x-axis of coordinate system
394 const vector e0 =
396 (
397 n ^ (points[prevE.otherVertex(prevVerti)] - points[prevVerti])
398 );
399
400 // y-axis of coordinate system
401 const vector e1 = n ^ e0;
402
403 if (mag(mag(e1) - 1) > SMALL)
404 {
405 {
406 Pout<< "Writing face:" << facei << " to face.obj" << endl;
407 OFstream str("face.obj");
408 writeOBJ(points, edges, fEdges, str);
409
410 Pout<< "Writing connectedEdges edges to faceEdges.obj" << endl;
411 writeLocalOBJ(points, edges, connectedEdges, "faceEdges.obj");
412 }
413
415 << "Unnormalized normal e1:" << e1
416 << " formed from cross product of e0:" << e0 << " n:" << n
417 << abort(FatalError);
418 }
419
420
421 //
422 // Determine maximum angle over all connected (unvisited) edges.
423 //
424
425 scalar maxAngle = -GREAT;
426 label maxEdgei = -1;
427
428 for (const label edgei : connectedEdges)
429 {
430 if (edgei != prevEdgei)
431 {
432 label stat = visited[edgei];
433
434 const edge& e = edges[edgei];
435
436 // Find out whether walk of edge from prevVert would be acceptable.
437 if
438 (
439 stat == UNVISITED
440 || (
441 stat == STARTTOEND
442 && prevVerti == e[1]
443 )
444 || (
445 stat == ENDTOSTART
446 && prevVerti == e[0]
447 )
448 )
449 {
450 // Calculate angle of edge with respect to base e0, e1
451 const vector vec =
453 (
454 n
455 ^ (points[e.otherVertex(prevVerti)] - points[prevVerti])
456 );
457
458 scalar angle = pseudoAngle(e0, e1, vec);
459
460 if (angle > maxAngle)
461 {
462 maxAngle = angle;
463 maxEdgei = edgei;
464 }
465 }
466 }
467 }
468
469
470 if (maxEdgei == -1)
471 {
472 // No unvisited edge found
473 {
474 Pout<< "Writing face:" << facei << " to face.obj" << endl;
475 OFstream str("face.obj");
476 writeOBJ(points, edges, fEdges, str);
477
478 Pout<< "Writing connectedEdges edges to faceEdges.obj" << endl;
479 writeLocalOBJ(points, edges, connectedEdges, "faceEdges.obj");
480 }
481
483 << "Trying to step from edge " << edges[prevEdgei]
484 << ", vertex " << prevVerti
485 << " but cannot find 'unvisited' edges among candidates:"
486 << connectedEdges
487 << abort(FatalError);
488 }
489
490 if (debug & 2)
491 {
492 Pout<< "Stepped from edge:" << edges[prevEdgei]
493 << " to edge:" << maxEdgei << " angle:" << edges[maxEdgei]
494 << " with angle:" << maxAngle
495 << endl;
496 }
497
498 return maxEdgei;
499}
500
501
502Foam::face Foam::intersectedSurface::walkFace
503(
504 const edgeSurface& eSurf,
505 const label facei,
506 const vector& n,
507 const Map<DynamicList<label>>& facePointEdges,
508
509 const label startEdgei,
510 const label startVerti,
511
512 Map<label>& visited
513)
514{
515 const pointField& points = eSurf.points();
516 const edgeList& edges = eSurf.edges();
517
518 // Overestimate size of face
519 face f(eSurf.faceEdges()[facei].size());
520
521 label fp = 0;
522
523 label verti = startVerti;
524 label edgei = startEdgei;
525
526 while (true)
527 {
528 const edge& e = edges[edgei];
529
530 if (debug & 2)
531 {
532 Pout<< "Now at:" << endl
533 << " edge:" << edgei << " vertices:" << e
534 << " positions:" << points[e.start()] << ' ' << points[e.end()]
535 << " vertex:" << verti << endl;
536 }
537
538 // Mark edge as visited
539 if (e[0] == verti)
540 {
541 visited[edgei] |= STARTTOEND;
542 }
543 else
544 {
545 visited[edgei] |= ENDTOSTART;
546 }
547
548
549 // Store face vertex
550 f[fp++] = verti;
551
552 // Step to other vertex
553 verti = e.otherVertex(verti);
554
555 if (verti == startVerti)
556 {
557 break;
558 }
559
560 // Step from vertex to next edge
561 edgei = nextEdge
562 (
563 eSurf,
564 visited,
565 facei,
566 n,
567 facePointEdges,
568 edgei,
569 verti
570 );
571 }
572
573 f.setSize(fp);
574
575 return f;
576}
577
578
579void Foam::intersectedSurface::findNearestVisited
580(
581 const edgeSurface& eSurf,
582 const label facei,
583 const Map<DynamicList<label>>& facePointEdges,
584 const Map<label>& pointVisited,
585 const point& pt,
586 const label excludePointi,
587
588 label& minVerti,
589 scalar& minDist
590)
591{
592 minVerti = -1;
593 minDist = GREAT;
594
595 forAllConstIters(pointVisited, iter)
596 {
597 const label pointi = iter.key();
598 const label nVisits = iter.val();
599
600 if (pointi != excludePointi)
601 {
602 if (nVisits == 2*facePointEdges[pointi].size())
603 {
604 // Fully visited (i.e. both sides of all edges)
605 scalar dist = mag(eSurf.points()[pointi] - pt);
606
607 if (dist < minDist)
608 {
609 minDist = dist;
610 minVerti = pointi;
611 }
612 }
613 }
614 }
615
616 if (minVerti == -1)
617 {
618 const labelList& fEdges = eSurf.faceEdges()[facei];
619
621 << "Dumping face edges to faceEdges.obj" << endl;
622
623 writeLocalOBJ(eSurf.points(), eSurf.edges(), fEdges, "faceEdges.obj");
624
626 << "No fully visited edge found for pt " << pt
627 << abort(FatalError);
628 }
629}
630
631
632Foam::faceList Foam::intersectedSurface::resplitFace
633(
634 const triSurface& surf,
635 const label facei,
636 const Map<DynamicList<label>>& facePointEdges,
637 const Map<label>& visited,
638 edgeSurface& eSurf
639)
640{
641 // Count the number of times point has been visited so we
642 // can compare number to facePointEdges.
643 Map<label> pointVisited(2*facePointEdges.size());
644
645 forAllConstIters(visited, iter)
646 {
647 const label edgei = iter.key();
648 const label stat = iter.val();
649
650 const edge& e = eSurf.edges()[edgei];
651
652 if (stat == STARTTOEND || stat == ENDTOSTART)
653 {
654 incCount(pointVisited, e[0], 1);
655 incCount(pointVisited, e[1], 1);
656 }
657 else if (stat == BOTH)
658 {
659 incCount(pointVisited, e[0], 2);
660 incCount(pointVisited, e[1], 2);
661 }
662 else if (stat == UNVISITED)
663 {
664 incCount(pointVisited, e[0], 0);
665 incCount(pointVisited, e[1], 0);
666 }
667 }
668
669
670 if (debug)
671 {
672 forAllConstIters(pointVisited, iter)
673 {
674 const label pointi = iter.key();
675 const label nVisits = iter.val();
676
677 Pout<< "point:" << pointi
678 << " nVisited:" << nVisits
679 << " pointEdges:" << facePointEdges[pointi].size() << endl;
680 }
681 }
682
683
684 // Find nearest point pair where one is not fully visited and
685 // the other is.
686
687 label visitedVert0 = -1;
688 label unvisitedVert0 = -1;
689
690 {
691 scalar minDist = GREAT;
692
693 forAllConstIters(facePointEdges, iter)
694 {
695 const label pointi = iter.key();
696
697 const DynamicList<label>& pEdges = iter.val();
698
699 const label nVisits = pointVisited[pointi];
700
701 if (nVisits < 2*pEdges.size())
702 {
703 // Not fully visited. Find nearest fully visited.
704
705 scalar nearDist;
706 label nearVerti;
707
708 findNearestVisited
709 (
710 eSurf,
711 facei,
712 facePointEdges,
713 pointVisited,
714 eSurf.points()[pointi],
715 -1, // Do not exclude vertex
716 nearVerti,
717 nearDist
718 );
719
720
721 if (nearDist < minDist)
722 {
723 minDist = nearDist;
724 visitedVert0 = nearVerti;
725 unvisitedVert0 = pointi;
726 }
727 }
728 }
729 }
730
731
732 // Find second intersection.
733 {
734 scalar minDist = GREAT;
735
736 forAllConstIters(facePointEdges, iter)
737 {
738 const label pointi = iter.key();
739
740 const DynamicList<label>& pEdges = iter.val();
741
742 if (pointi != unvisitedVert0)
743 {
744 const label nVisits = pointVisited[pointi];
745
746 if (nVisits < 2*pEdges.size())
747 {
748 // Not fully visited. Find nearest fully visited.
749
750 scalar nearDist;
751 label nearVerti;
752
753 findNearestVisited
754 (
755 eSurf,
756 facei,
757 facePointEdges,
758 pointVisited,
759 eSurf.points()[pointi],
760 visitedVert0, // vertex to exclude
761 nearVerti,
762 nearDist
763 );
764
765
766 if (nearDist < minDist)
767 {
768 minDist = nearDist;
769 }
770 }
771 }
772 }
773 }
774
775
776 // Add the new intersection edges to the edgeSurface
777 edgeList additionalEdges(1);
778 additionalEdges[0] = edge(visitedVert0, unvisitedVert0);
779
780 eSurf.addIntersectionEdges(facei, additionalEdges);
781
782 if (debug)
783 {
784 fileName newFName("face_" + Foam::name(facei) + "_newEdges.obj");
785 Pout<< "Dumping face:" << facei << " to " << newFName << endl;
786 writeLocalOBJ
787 (
788 eSurf.points(),
789 eSurf.edges(),
790 eSurf.faceEdges()[facei],
791 newFName
792 );
793 }
794
795 // Retry splitFace. Use recursion since is rare situation.
796 return splitFace(surf, facei, eSurf);
797}
798
799
800Foam::faceList Foam::intersectedSurface::splitFace
801(
802 const triSurface& surf,
803 const label facei,
804 edgeSurface& eSurf
805)
806{
807 // Alias
808 const pointField& points = eSurf.points();
809 const edgeList& edges = eSurf.edges();
810 const labelList& fEdges = eSurf.faceEdges()[facei];
811
812 // Create local (for the face only) point-edge connectivity.
813 Map<DynamicList<label>> facePointEdges
814 (
815 calcPointEdgeAddressing
816 (
817 eSurf,
818 facei
819 )
820 );
821
822 // Order in which edges have been walked. Initialize outside edges.
823 Map<label> visited(fEdges.size()*2);
824
825 forAll(fEdges, i)
826 {
827 label edgei = fEdges[i];
828
829 if (eSurf.isSurfaceEdge(edgei))
830 {
831 // Edge is edge from original surface so an outside edge for
832 // the current face.
833 label surfEdgei = eSurf.parentEdge(edgei);
834
835 label owner = surf.edgeOwner()[surfEdgei];
836
837 if
838 (
839 owner == facei
840 || sameEdgeOrder
841 (
842 surf.localFaces()[owner],
843 surf.localFaces()[facei]
844 )
845 )
846 {
847 // Edge is in same order as current face.
848 // Mark off the opposite order.
849 visited.insert(edgei, ENDTOSTART);
850 }
851 else
852 {
853 // Edge is in same order as current face.
854 // Mark off the opposite order.
855 visited.insert(edgei, STARTTOEND);
856 }
857 }
858 else
859 {
860 visited.insert(edgei, UNVISITED);
861 }
862 }
863
864
865
866 // Storage for faces
867 DynamicList<face> faces(fEdges.size());
868
869 while (true)
870 {
871 // Find starting edge:
872 // - unvisited triangle edge
873 // - once visited intersection edge
874 // Give priority to triangle edges.
875 label startEdgei = -1;
876 label startVerti = -1;
877
878 forAll(fEdges, i)
879 {
880 label edgei = fEdges[i];
881
882 const edge& e = edges[edgei];
883
884 label stat = visited[edgei];
885
886 if (stat == STARTTOEND)
887 {
888 startEdgei = edgei;
889 startVerti = e[1];
890
891 if (eSurf.isSurfaceEdge(edgei))
892 {
893 break;
894 }
895 }
896 else if (stat == ENDTOSTART)
897 {
898 startEdgei = edgei;
899 startVerti = e[0];
900
901 if (eSurf.isSurfaceEdge(edgei))
902 {
903 break;
904 }
905 }
906 }
907
908 if (startEdgei == -1)
909 {
910 break;
911 }
912
913 //Pout<< "splitFace: starting walk from edge:" << startEdgei
914 // << ' ' << edges[startEdgei] << " vertex:" << startVerti << endl;
915
917 //printVisit(eSurf.edges(), fEdges, visited);
918
919 //{
920 // Pout<< "Writing face:" << faceI << " to face.obj" << endl;
921 // OFstream str("face.obj");
922 // writeOBJ(eSurf.points(), eSurf.edges(), fEdges, str);
923 //}
924
925 faces.append
926 (
927 walkFace
928 (
929 eSurf,
930 facei,
931 surf.faceNormals()[facei],
932 facePointEdges,
933
934 startEdgei,
935 startVerti,
936
937 visited
938 )
939 );
940 }
941
942
943 // Check if any unvisited edges left.
944 forAll(fEdges, i)
945 {
946 label edgei = fEdges[i];
947
948 label stat = visited[edgei];
949
950 if (eSurf.isSurfaceEdge(edgei) && stat != BOTH)
951 {
953 << "Dumping face edges to faceEdges.obj" << endl;
954
955 writeLocalOBJ(points, edges, fEdges, "faceEdges.obj");
956
958 << "Problem: edge " << edgei << " vertices "
959 << edges[edgei] << " on face " << facei
960 << " has visited status " << stat << " from a "
961 << "right-handed walk along all"
962 << " of the triangle edges. Are the original surfaces"
963 << " closed and non-intersecting?"
964 << abort(FatalError);
965 }
966 else if (stat != BOTH)
967 {
968 // Redo face after introducing extra edge. Edge introduced
969 // should be one nearest to any fully visited edge.
970 return resplitFace
971 (
972 surf,
973 facei,
974 facePointEdges,
975 visited,
976 eSurf
977 );
978 }
979 }
980
981
982 // See if normal needs flipping.
983 faces.shrink();
984
985 const vector n = faces[0].areaNormal(eSurf.points());
986
987 if ((n & surf.faceNormals()[facei]) < 0)
988 {
989 for (face& f : faces)
990 {
991 reverse(f);
992 }
993 }
994
995 return faces;
996}
997
998
999// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1000
1002:
1003 triSurface(),
1004 intersectionEdges_(0),
1005 faceMap_(0),
1006 nSurfacePoints_(0)
1007{}
1008
1009
1011:
1012 triSurface(surf),
1013 intersectionEdges_(0),
1014 faceMap_(0),
1015 nSurfacePoints_(0)
1016{}
1017
1018
1020(
1021 const triSurface& surf,
1022 const bool isFirstSurface,
1023 const surfaceIntersection& inter
1024)
1025:
1026 triSurface(),
1027 intersectionEdges_(0),
1028 faceMap_(0),
1029 nSurfacePoints_(surf.nPoints())
1030{
1031 if (inter.cutPoints().empty() && inter.cutEdges().empty())
1032 {
1033 // No intersection. Make straight copy.
1035
1036 // Identity for face map
1037 faceMap_.setSize(size());
1038
1039 forAll(faceMap_, facei)
1040 {
1041 faceMap_[facei] = facei;
1042 }
1043 return;
1044 }
1045
1046
1047 // Calculate face-edge addressing for surface + intersection.
1048 edgeSurface eSurf(surf, isFirstSurface, inter);
1049
1050 // Update point information for any removed points. (when are they removed?
1051 // - but better make sure)
1052 nSurfacePoints_ = eSurf.nSurfacePoints();
1053
1054 // Now we have full points, edges and edge addressing for surf. Start
1055 // extracting faces and triangulate them.
1056
1057 // Storage for new triangles of surface.
1058 DynamicList<labelledTri> newTris(eSurf.edges().size()/2);
1059
1060 // Start in newTris for decomposed face.
1061 labelList startTrii(surf.size(), Zero);
1062
1063 forAll(surf, facei)
1064 {
1065 startTrii[facei] = newTris.size();
1066
1067 if (eSurf.faceEdges()[facei].size() != surf.faceEdges()[facei].size())
1068 {
1069 // Face has been cut by intersection.
1070 // Cut face into multiple subfaces. Use faceEdge information
1071 // from edgeSurface only. (original triSurface 'surf' is used for
1072 // faceNormals and region number only)
1073 faceList newFaces
1074 (
1075 splitFace
1076 (
1077 surf,
1078 facei, // current triangle
1079 eSurf // face-edge description of surface
1080 // + intersection
1081 )
1082 );
1083 forAll(newFaces, newFacei)
1084 {
1085 const face& newF = newFaces[newFacei];
1086 const vector& n = surf.faceNormals()[facei];
1087 const label region = surf[facei].region();
1088
1089 faceTriangulation tris(eSurf.points(), newF, n);
1090
1091 forAll(tris, trii)
1092 {
1093 const triFace& t = tris[trii];
1094
1095 forAll(t, i)
1096 {
1097 if (t[i] < 0 || t[i] >= eSurf.points().size())
1098 {
1100 << "Face triangulation of face " << facei
1101 << " uses points outside range 0.."
1102 << eSurf.points().size()-1 << endl
1103 << "Triangulation:"
1104 << tris << abort(FatalError);
1105 }
1106 }
1107
1108 newTris.append(labelledTri(t[0], t[1], t[2], region));
1109 }
1110 }
1111 }
1112 else
1113 {
1114 // Face has not been cut at all. No need to renumber vertices since
1115 // eSurf keeps surface vertices first.
1116 newTris.append(surf.localFaces()[facei]);
1117 }
1118 }
1119
1120 newTris.shrink();
1121
1122
1123 // Construct triSurface. Note that addressing will be compact since
1124 // edgeSurface is compact.
1126 (
1128 (
1129 newTris,
1130 surf.patches(),
1131 eSurf.points()
1132 )
1133 );
1134
1135 // Construct mapping back into original surface
1136 faceMap_.setSize(size());
1137
1138 for (label facei = 0; facei < surf.size()-1; facei++)
1139 {
1140 for (label trii = startTrii[facei]; trii < startTrii[facei+1]; trii++)
1141 {
1142 faceMap_[trii] = facei;
1143 }
1144 }
1145 for (label trii = startTrii[surf.size()-1]; trii < size(); trii++)
1146 {
1147 faceMap_[trii] = surf.size()-1;
1148 }
1149
1150
1151 // Find edges on *this which originate from 'cuts'. (i.e. newEdgei >=
1152 // nSurfaceEdges) Renumber edges into local triSurface numbering.
1153
1154 intersectionEdges_.setSize(eSurf.edges().size() - eSurf.nSurfaceEdges());
1155
1156 label intersectionEdgei = 0;
1157
1158 for
1159 (
1160 label edgei = eSurf.nSurfaceEdges();
1161 edgei < eSurf.edges().size();
1162 edgei++
1163 )
1164 {
1165 // Get edge vertices in triSurface local numbering
1166 const edge& e = eSurf.edges()[edgei];
1167 label surfStarti = meshPointMap()[e.start()];
1168 label surfEndi = meshPointMap()[e.end()];
1169
1170 // Find edge connected to surfStarti which also uses surfEndi.
1171 label surfEdgei = -1;
1172
1173 const labelList& pEdges = pointEdges()[surfStarti];
1174
1175 forAll(pEdges, i)
1176 {
1177 const edge& surfE = edges()[pEdges[i]];
1178
1179 // Edge already connected to surfStart for sure. See if also
1180 // connects to surfEnd
1181 if (surfE.found(surfEndi))
1182 {
1183 surfEdgei = pEdges[i];
1184
1185 break;
1186 }
1187 }
1188
1189 if (surfEdgei != -1)
1190 {
1191 intersectionEdges_[intersectionEdgei++] = surfEdgei;
1192 }
1193 else
1194 {
1196 << "Cannot find edge among candidates " << pEdges
1197 << " which uses points " << surfStarti
1198 << " and " << surfEndi
1199 << abort(FatalError);
1200 }
1201 }
1202}
1203
1204
1205// ************************************************************************* //
label n
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:434
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
const labelListList & pointEdges() const
Return point-edge addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const Map< label > & meshPointMap() const
Mesh point map.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
const labelListList & faceEdges() const
Return face-edge addressing.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
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 size() const noexcept
The number of elements in the UList.
Definition: UListI.H:420
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Description of surface in form of 'cloud of edges'.
Definition: edgeSurface.H:76
label nSurfaceEdges() const
Definition: edgeSurface.H:141
const edgeList & edges() const
Definition: edgeSurface.H:136
label nSurfacePoints() const
Definition: edgeSurface.H:131
const pointField & points() const
Definition: edgeSurface.H:126
const labelListList & faceEdges() const
From face to our edges_.
Definition: edgeSurface.H:170
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
bool found(const label pointLabel) const
Return true if point label is found in edge.
Definition: edgeI.H:136
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
Triangulation of faces. Handles concave polygons as well (inefficiently)
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Given triSurface and intersection creates the intersected (properly triangulated) surface....
intersectedSurface()
Construct null.
static const label ENDTOSTART
static const label UNVISITED
static const label STARTTOEND
A triFace with additional (region) index.
Definition: labelledTri.H:60
Basic surface-surface intersection description. Constructed from two surfaces it creates a descriptio...
const edgeList & cutEdges() const
The list of created edges.
const pointField & cutPoints() const
The list of cut points.
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:72
Triangulated surface description with patch information.
Definition: triSurface.H:79
triSurface()
Default construct.
Definition: triSurface.C:432
const geometricSurfacePatchList & patches() const noexcept
Definition: triSurface.H:399
void operator=(const triSurface &surf)
Copy assignment.
Definition: triSurface.C:999
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
label nPoints
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Namespace for OpenFOAM.
scalar pseudoAngle(const vector &e0, const vector &e1, const vector &vec)
Estimate angle of vec in coordinate system (e0, e1, e0^e1).
Definition: transform.H:401
List< label > labelList
A List of labels.
Definition: List.H:66
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition: UListI.H:449
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
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.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:260
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
3D tensor transformation operations.