Loading...
Searching...
No Matches
surfaceBooleanFeatures.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) 2016-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Application
28 surfaceBooleanFeatures
29
30Group
31 grpSurfaceUtilities
32
33Description
34 Generates the extendedFeatureEdgeMesh for the interface between a boolean
35 operation on two surfaces.
36
37 Assumes that the orientation of the surfaces is correct:
38 - if the operation is union or intersection, that both surface's normals
39 (n) have the same orientation with respect to a point, i.e. surfaces
40 A and B are orientated the same with respect to point x:
41
42 \verbatim
43 _______
44 | |--> n
45 | ___|___ x
46 |A | | |--> n
47 |___|___| B|
48 | |
49 |_______|
50
51 \endverbatim
52
53 - if the operation is a subtraction, the surfaces should be oppositely
54 oriented with respect to a point, i.e. for (A - B), then B's orientation
55 should be such that x is "inside", and A's orientation such that x is
56 "outside"
57
58 \verbatim
59 _______
60 | |--> n
61 | ___|___ x
62 |A | | |
63 |___|___| B|
64 | n <--|
65 |_______|
66
67 \endverbatim
68
69 When the operation is performed - for union, all of the edges generated
70 where one surfaces cuts another are all "internal" for union,
71 and "external" for intersection, (B - A) and (A - B).
72 This has been assumed, formal (dis)proof is invited.
73
74\*---------------------------------------------------------------------------*/
75
76#include "triSurface.H"
77#include "argList.H"
78#include "Time.H"
79#include "featureEdgeMesh.H"
81#include "triSurfaceSearch.H"
82#include "triSurfaceMesh.H"
83#include "OFstream.H"
84#include "booleanSurface.H"
85#include "edgeIntersections.H"
86#include "meshTools.H"
87#include "DynamicField.H"
88#include "Enum.H"
89
90#ifndef NO_CGAL
91
92// Silence boost bind deprecation warnings (before CGAL-5.2.1)
93#include "CGAL/version.h"
94#if defined(CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1050211000)
95#define BOOST_BIND_GLOBAL_PLACEHOLDERS
96#endif
97#pragma clang diagnostic ignored "-Wbitwise-instead-of-logical"
98#pragma clang diagnostic ignored "-Wdeprecated-builtins"
99#pragma clang diagnostic ignored "-Wdeprecated-declarations"
100
101#include <CGAL/AABB_tree.h>
102#if defined(CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1060011000)
103#include <CGAL/AABB_traits.h>
104#else
105#include <CGAL/AABB_traits_3.h>
106#endif
107#include <CGAL/AABB_face_graph_triangle_primitive.h>
109#include "PolyhedronReader.H"
110typedef CGAL::AABB_face_graph_triangle_primitive
111<
112 Polyhedron, CGAL::Default, CGAL::Tag_false
113> Primitive;
114#if defined(CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1060011000)
115typedef CGAL::AABB_traits<K, Primitive> Traits;
116#else
117typedef CGAL::AABB_traits_3<K, Primitive> Traits;
118#endif
119typedef CGAL::AABB_tree<Traits> Tree;
120
121// Used boost::optional prior to CGAL-6.0
122#if (CGAL_VERSION_NR < 1060000910)
123typedef boost::optional
124#else
125typedef std::optional
126#endif
127<
128 Tree::Intersection_and_primitive_id<Segment>::Type
129> Segment_intersection;
130
131#endif // NO_CGAL
132
133
134using namespace Foam;
135
136// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137
138bool intersectSurfaces
139(
140 triSurface& surf,
141 edgeIntersections& edgeCuts
142)
143{
144 bool hasMoved = false;
145
146 for (label iter = 0; iter < 10; iter++)
147 {
148 Info<< "Determining intersections of surface edges with itself" << endl;
149
150 // Determine surface edge intersections. Allow surface to be moved.
151
152 // Number of iterations needed to resolve degenerates
153 label nIters = 0;
154 {
155 triSurfaceSearch querySurf(surf);
156
157 scalarField surfPointTol
158 (
159 max(1e-3*edgeIntersections::minEdgeLength(surf), SMALL)
160 );
161
162 // Determine raw intersections
163 edgeCuts = edgeIntersections
164 (
165 surf,
166 querySurf,
167 surfPointTol
168 );
169
170 // Shuffle a bit to resolve degenerate edge-face hits
171 {
172 pointField points(surf.points());
173
174 nIters =
175 edgeCuts.removeDegenerates
176 (
177 5, // max iterations
178 surf,
179 querySurf,
180 surfPointTol,
181 points // work array
182 );
183
184 if (nIters != 0)
185 {
186 // Update geometric quantities
187 surf.movePoints(points);
188 hasMoved = true;
189 }
190 }
191 }
192 }
193
194 if (hasMoved)
195 {
196 fileName newFile("surf.obj");
197 Info<< "Surface has been moved. Writing to " << newFile << endl;
198 surf.write(newFile);
199 }
200
201 return hasMoved;
202}
203
204
205// Keep on shuffling surface points until no more degenerate intersections.
206// Moves both surfaces and updates set of edge cuts.
207bool intersectSurfaces
208(
209 triSurface& surf1,
210 edgeIntersections& edgeCuts1,
211 triSurface& surf2,
212 edgeIntersections& edgeCuts2
213)
214{
215 bool hasMoved1 = false;
216 bool hasMoved2 = false;
217
218 for (label iter = 0; iter < 10; iter++)
219 {
220 Info<< "Determining intersections of surf1 edges with surf2"
221 << " faces" << endl;
222
223 // Determine surface1 edge intersections. Allow surface to be moved.
224
225 // Number of iterations needed to resolve degenerates
226 label nIters1 = 0;
227 {
228 triSurfaceSearch querySurf2(surf2);
229
230 scalarField surf1PointTol
231 (
232 max(1e-3*edgeIntersections::minEdgeLength(surf1), SMALL)
233 );
234
235 // Determine raw intersections
236 edgeCuts1 = edgeIntersections
237 (
238 surf1,
239 querySurf2,
240 surf1PointTol
241 );
242
243 // Shuffle a bit to resolve degenerate edge-face hits
244 {
245 pointField points1(surf1.points());
246
247 nIters1 =
248 edgeCuts1.removeDegenerates
249 (
250 5, // max iterations
251 surf1,
252 querySurf2,
253 surf1PointTol,
254 points1 // work array
255 );
256
257 if (nIters1 != 0)
258 {
259 // Update geometric quantities
260 surf1.movePoints(points1);
261 hasMoved1 = true;
262 }
263 }
264 }
265
266 Info<< "Determining intersections of surf2 edges with surf1"
267 << " faces" << endl;
268
269 label nIters2 = 0;
270 {
271 triSurfaceSearch querySurf1(surf1);
272
273 scalarField surf2PointTol
274 (
275 max(1e-3*edgeIntersections::minEdgeLength(surf2), SMALL)
276 );
277
278 // Determine raw intersections
279 edgeCuts2 = edgeIntersections
280 (
281 surf2,
282 querySurf1,
283 surf2PointTol
284 );
285
286 // Shuffle a bit to resolve degenerate edge-face hits
287 {
288 pointField points2(surf2.points());
289
290 nIters2 =
291 edgeCuts2.removeDegenerates
292 (
293 5, // max iterations
294 surf2,
295 querySurf1,
296 surf2PointTol,
297 points2 // work array
298 );
299
300 if (nIters2 != 0)
301 {
302 // Update geometric quantities
303 surf2.movePoints(points2);
304 hasMoved2 = true;
305 }
306 }
307 }
308
309 if (nIters1 == 0 && nIters2 == 0)
310 {
311 Info<< "** Resolved all intersections to be proper edge-face pierce"
312 << endl;
313 break;
314 }
315 }
316
317 if (hasMoved1)
318 {
319 fileName newFile("surf1.obj");
320 Info<< "Surface 1 has been moved. Writing to " << newFile
321 << endl;
322 surf1.write(newFile);
323 }
324
325 if (hasMoved2)
326 {
327 fileName newFile("surf2.obj");
328 Info<< "Surface 2 has been moved. Writing to " << newFile
329 << endl;
330 surf2.write(newFile);
331 }
332
333 return hasMoved1 || hasMoved2;
334}
335
336
337label calcNormalDirection
338(
339 const vector& normal,
340 const vector& otherNormal,
341 const vector& edgeDir,
342 const vector& faceCentre,
343 const vector& pointOnEdge
344)
345{
346 const vector cross = normalised(normal ^ edgeDir);
347
348 const vector fC0tofE0 = normalised(faceCentre - pointOnEdge);
349
350 label nDir = ((cross & fC0tofE0) > 0.0 ? 1 : -1);
351
352 nDir *= ((otherNormal & fC0tofE0) > 0.0 ? -1 : 1);
353
354 return nDir;
355}
356
357
358void calcEdgeCuts
359(
360 triSurface& surf1,
361 triSurface& surf2,
362 const bool perturb,
363 edgeIntersections& edgeCuts1,
364 edgeIntersections& edgeCuts2
365)
366{
367 if (perturb)
368 {
369 intersectSurfaces
370 (
371 surf1,
372 edgeCuts1,
373 surf2,
374 edgeCuts2
375 );
376 }
377 else
378 {
379 triSurfaceSearch querySurf2(surf2);
380
381 Info<< "Determining intersections of surf1 edges with surf2 faces"
382 << endl;
383
384 edgeCuts1 = edgeIntersections
385 (
386 surf1,
387 querySurf2,
388 max(1e-3*edgeIntersections::minEdgeLength(surf1), SMALL)
389 );
390
391 triSurfaceSearch querySurf1(surf1);
392
393 Info<< "Determining intersections of surf2 edges with surf1 faces"
394 << endl;
395
396 edgeCuts2 = edgeIntersections
397 (
398 surf2,
399 querySurf1,
400 max(1e-3*edgeIntersections::minEdgeLength(surf2), SMALL)
401 );
402 }
403}
404
405
406// CGAL variants
407
408#ifndef NO_CGAL
409
410void visitPointRegion
411(
412 const triSurface& s,
413 const label zoneI,
414 const label pointI,
415 const label startEdgeI,
416 const label startFaceI,
417 labelList& pFacesZone
418)
419{
420 const labelList& eFaces = s.edgeFaces()[startEdgeI];
421
422 if (eFaces.size() == 2)
423 {
424 label nextFaceI;
425 if (eFaces[0] == startFaceI)
426 {
427 nextFaceI = eFaces[1];
428 }
429 else if (eFaces[1] == startFaceI)
430 {
431 nextFaceI = eFaces[0];
432 }
433 else
434 {
436 << "problem" << exit(FatalError);
437
438 nextFaceI = -1;
439 }
440
441
442
443 label index = s.pointFaces()[pointI].find(nextFaceI);
444
445 if (pFacesZone[index] == -1)
446 {
447 // Mark face as been visited.
448 pFacesZone[index] = zoneI;
449
450 // Step to next edge on face which is still using pointI
451 const labelList& fEdges = s.faceEdges()[nextFaceI];
452
453 label nextEdgeI = -1;
454
455 forAll(fEdges, i)
456 {
457 label edgeI = fEdges[i];
458 const edge& e = s.edges()[edgeI];
459
460 if (edgeI != startEdgeI && (e[0] == pointI || e[1] == pointI))
461 {
462 nextEdgeI = edgeI;
463
464 break;
465 }
466 }
467
468 if (nextEdgeI == -1)
469 {
471 << "Problem: cannot find edge out of " << fEdges
472 << "on face " << nextFaceI << " that uses point " << pointI
473 << " and is not edge " << startEdgeI << abort(FatalError);
474 }
475
476
477 visitPointRegion
478 (
479 s,
480 zoneI,
481 pointI,
482 nextEdgeI,
483 nextFaceI,
484 pFacesZone
485 );
486 }
487 }
488}
489
490
491label dupNonManifoldPoints(triSurface& s, labelList& pointMap)
492{
493 const labelListList& pf = s.pointFaces();
494 const labelListList& fe = s.faceEdges();
495 const edgeList& edges = s.edges();
496
497
498 DynamicField<point> newPoints(s.points());
499 // From dupSurf back to s.pointa
500 DynamicList<label> newPointMap(identity(newPoints.size()));
501 List<labelledTri> newFaces(s);
502 label nNonManifold = 0;
503
504 forAll(pf, pointI)
505 {
506 const labelList& pFaces = pf[pointI];
507
508 // Visited faces (as indices into pFaces)
509 labelList pFacesZone(pFaces.size(), -1);
510
511 label nZones = 0;
512 label index = 0;
513 for (; index < pFacesZone.size(); index++)
514 {
515 if (pFacesZone[index] == -1)
516 {
517 label zoneI = nZones++;
518 pFacesZone[index] = zoneI;
519
520 label faceI = pFaces[index];
521 const labelList& fEdges = fe[faceI];
522
523 // Find starting edge
524 forAll(fEdges, fEdgeI)
525 {
526 label edgeI = fEdges[fEdgeI];
527 const edge& e = edges[edgeI];
528
529 if (e[0] == pointI || e[1] == pointI)
530 {
531 visitPointRegion
532 (
533 s,
534 zoneI,
535 pointI,
536 edgeI,
537 faceI,
538 pFacesZone
539 );
540 }
541 }
542 }
543 }
544
545
546 // Subset
547 if (nZones > 1)
548 {
549 for (label zoneI = 1; zoneI < nZones; zoneI++)
550 {
551 label newPointI = newPoints.size();
552 newPointMap.append(s.meshPoints()[pointI]);
553 newPoints.append(s.points()[s.meshPoints()[pointI]]);
554
555 forAll(pFacesZone, index)
556 {
557 if (pFacesZone[index] == zoneI)
558 {
559 label faceI = pFaces[index];
560 const labelledTri& localF = s.localFaces()[faceI];
561 forAll(localF, fp)
562 {
563 if (localF[fp] == pointI)
564 {
565 newFaces[faceI][fp] = newPointI;
566 }
567 }
568 }
569 }
570 }
571 nNonManifold++;
572 }
573 }
574
575
576 Info<< "Duplicating " << nNonManifold << " points out of " << s.nPoints()
577 << endl;
578 if (nNonManifold > 0)
579 {
580 triSurface dupSurf(newFaces, s.patches(), newPoints, true);
581
582 // Create map from dupSurf localPoints to s.localPoints
583 const Map<label> mpm = s.meshPointMap();
584
585 const labelList& dupMp = dupSurf.meshPoints();
586
587 labelList dupPointMap(dupMp.size());
588 forAll(dupMp, pointI)
589 {
590 label dupMeshPointI = dupMp[pointI];
591 label meshPointI = newPointMap[dupMeshPointI];
592 dupPointMap[pointI] = mpm[meshPointI];
593 }
594
595
596 forAll(dupPointMap, pointI)
597 {
598 const point& dupPt = dupSurf.points()[dupMp[pointI]];
599 label sLocalPointI = dupPointMap[pointI];
600 label sMeshPointI = s.meshPoints()[sLocalPointI];
601 const point& sPt = s.points()[sMeshPointI];
602
603 if (mag(dupPt-sPt) > SMALL)
604 {
606 << "dupPt:" << dupPt
607 << " sPt:" << sPt
608 << exit(FatalError);
609 }
610 }
611
612
613 //s.transfer(dupSurf);
614 s = dupSurf;
615 pointMap = labelUIndList(pointMap, dupPointMap)();
616 }
617
618 return nNonManifold;
619}
620
621
622// Find intersections of surf1 by edges of surf2. Return number of degenerate
623// cuts (cuts through faces/edges/points)
624labelPair edgeIntersectionsCGAL
625(
626 const Tree& tree,
627 const triSurface& surf,
628 const pointField& points,
629 edgeIntersections& edgeCuts
630)
631{
632 const edgeList& edges = surf.edges();
633 const labelList& meshPoints = surf.meshPoints();
634
635 //Info<< "Intersecting CGAL surface ..." << endl;
636 List<List<pointIndexHit>> intersections(edges.size());
637 labelListList classifications(edges.size());
638
639 label nPoints = 0;
640 label nSegments = 0;
641
642 std::vector<Segment_intersection> segments;
643 forAll(edges, edgeI)
644 {
645 const edge& e = edges[edgeI];
646
647 const point& a = points[meshPoints[e[0]]];
648 const point& b = points[meshPoints[e[1]]];
649
650 K::Segment_3 segment_query
651 (
652 Point(a.x(), a.y(), a.z()),
653 Point(b.x(), b.y(), b.z())
654 );
655
656 segments.clear();
657 tree.all_intersections(segment_query, std::back_inserter(segments));
658
659
660 for (const Segment_intersection& intersect : segments)
661 {
662 // Get intersection object
663 if
664 (
665 // Used boost::variant prior to CGAL-6.0
666 #if (CGAL_VERSION_NR < 1060000910)
667 const Point* ptPtr = boost::get<Point>(&(intersect->first))
668 #else
669 const Point* ptPtr = std::get_if<Point>(&(intersect->first))
670 #endif
671 )
672 {
673 point pt
674 (
675 CGAL::to_double(ptPtr->x()),
676 CGAL::to_double(ptPtr->y()),
677 CGAL::to_double(ptPtr->z())
678 );
679
680 #if defined (CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1041400000)
681 Polyhedron::Face_handle f = (intersect->second);
682 #else
683 // 1.14 and later
684 Polyhedron::Face_handle f = (intersect->second).first;
685 #endif
686
687 intersections[edgeI].append
688 (
690 (
691 true,
692 pt,
693 f->index
694 )
695 );
696 // Intersection on edge interior
697 classifications[edgeI].append(-1);
698 ++nPoints;
699 }
700 else if
701 (
702 // Used boost::variant prior to CGAL-6.0
703 #if (CGAL_VERSION_NR < 1060000910)
704 const Segment* sPtr = boost::get<Segment>(&(intersect->first))
705 #else
706 const Segment* sPtr = std::get_if<Segment>(&(intersect->first))
707 #endif
708 )
709 {
710 #if defined (CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1041400000)
711 Polyhedron::Face_handle f = (intersect->second);
712 #else
713 // 1.14 and later
714 Polyhedron::Face_handle f = (intersect->second).first;
715 #endif
716
717 //std::cout
718 // << "intersection object is a segment:" << sPtr->source()
719 // << " " << sPtr->target() << std::endl;
720
721 //std::cout<< "triangle:" << f->index
722 // << " region:" << f->region << std::endl;
723
724 const point source
725 (
726 CGAL::to_double(sPtr->source().x()),
727 CGAL::to_double(sPtr->source().y()),
728 CGAL::to_double(sPtr->source().z())
729 );
730
731 const point target
732 (
733 CGAL::to_double(sPtr->target().x()),
734 CGAL::to_double(sPtr->target().y()),
735 CGAL::to_double(sPtr->target().z())
736 );
737
738 // Make up some intersection point
739 intersections[edgeI].append
740 (
742 (
743 true,
744 0.5*(source+target),
745 f->index
746 )
747 );
748 // Intersection aligned with face. Tbd: enums
749 classifications[edgeI].append(2);
750 ++nSegments;
751 }
752 }
753 }
754
755 edgeCuts = edgeIntersections(intersections, classifications);
756
757 return labelPair(nPoints, nSegments);
758}
759
760
761// Intersect edges of surf1 until no more degenerate intersections. Return
762// number of degenerates
763labelPair edgeIntersectionsAndShuffleCGAL
764(
765 Random& rndGen,
766 const triSurface& surf2,
767 const scalarField& surf1PointTol,
768 triSurface& surf1,
769 edgeIntersections& edgeCuts1
770)
771{
772 //Info<< "Constructing CGAL surface ..." << endl;
774 PolyhedronReader(surf2, p);
775
776
777 //Info<< "Constructing CGAL tree ..." << endl;
778 const Tree tree(p.facets_begin(), p.facets_end(), p);
779
780
781 labelPair cuts(0, 0);
782
783 // Update surface1 points until no longer intersecting
784 pointField surf1Points(surf1.points());
785
786 for (label iter = 0; iter < 5; iter++)
787 {
788 // See which edges of 1 intersect 2
789 Info<< "Determining intersections of " << surf1.nEdges()
790 << " edges with surface of " << label(tree.size()) << " triangles"
791 << endl;
792 cuts = edgeIntersectionsCGAL
793 (
794 tree,
795 surf1,
796 surf1Points,
797 edgeCuts1
798 );
799 Info<< "Determined intersections:" << nl
800 << " edges : " << surf1.nEdges() << nl
801 << " non-degenerate cuts : " << cuts.first() << nl
802 << " degenerate cuts : " << cuts.second() << nl
803 << endl;
804
805 if (cuts.second() == 0)
806 {
807 break;
808 }
809
810 Info<< "Shuffling conflicting points" << endl;
811
812 const labelListList& edgeStat = edgeCuts1.classification();
813 const edgeList& edges = surf1.edges();
814 const labelList& mp = surf1.meshPoints();
815 const point p05(0.5, 0.5, 0.5);
816
817 forAll(edgeStat, edgeI)
818 {
819 const labelList& stat = edgeStat[edgeI];
820 forAll(stat, i)
821 {
822 if (stat[i] == 2)
823 {
824 const edge& e = edges[edgeI];
825 forAll(e, eI)
826 {
827 vector d = rndGen.sample01<vector>() - p05;
828 surf1Points[mp[e[eI]]] += surf1PointTol[e[eI]]*d;
829 }
830 }
831 }
832 }
833 }
834 surf1.movePoints(surf1Points);
835
836 return cuts;
837}
838
839
840// Return map from subSurf.edges() back to surf.edges()
841labelList matchEdges
842(
843 const triSurface& surf,
844 const triSurface& subSurf,
845 const labelList& pointMap
846)
847{
848 if (pointMap.size() != subSurf.nPoints())
849 {
851 << "problem" << exit(FatalError);
852 }
853
854 labelList edgeMap(subSurf.nEdges(), -1);
855
856 const edgeList& edges = surf.edges();
857 const labelListList& pointEdges = surf.pointEdges();
858
859 const edgeList& subEdges = subSurf.edges();
860
861
862 forAll(subEdges, subEdgeI)
863 {
864 const edge& subE = subEdges[subEdgeI];
865
866 // Match points on edge to those on surf
867 label start = pointMap[subE[0]];
868 label end = pointMap[subE[1]];
869 const labelList& pEdges = pointEdges[start];
870 forAll(pEdges, pEdgeI)
871 {
872 label edgeI = pEdges[pEdgeI];
873 const edge& e = edges[edgeI];
874
875 if (e.otherVertex(start) == end)
876 {
877 if (edgeMap[subEdgeI] == -1)
878 {
879 edgeMap[subEdgeI] = edgeI;
880 }
881 else if (edgeMap[subEdgeI] != edgeI)
882 {
884 << subE << " points:"
885 << subE.line(subSurf.localPoints())
886 << " matches to " << edgeI
887 << " points:" << e.line(surf.localPoints())
888 << " but also " << edgeMap[subEdgeI]
889 << " points:"
890 << edges[edgeMap[subEdgeI]].line(surf.localPoints())
891 << exit(FatalError);
892 }
893 break;
894 }
895 }
896
897 if (edgeMap[subEdgeI] == -1)
898 {
900 << subE << " at:" << subSurf.localPoints()[subE[0]]
901 << subSurf.localPoints()[subE[1]]
902 << exit(FatalError);
903 }
904 }
905
906 return edgeMap;
907}
908
909
910void calcEdgeCutsCGAL
911(
912 triSurface& surf1,
913 triSurface& surf2,
914 const bool perturb,
915 edgeIntersections& edgeCuts1,
916 edgeIntersections& edgeCuts2
917)
918{
919 if (!perturb)
920 {
921 // See which edges of 1 intersect 2
922 {
923 Info<< "Intersect surface 1 edges with surface 2:" << nl;
924 Info<< " constructing CGAL surface ..." << endl;
926 PolyhedronReader(surf2, p);
927
928 Info<< " constructing CGAL tree ..." << endl;
929 const Tree tree(p.facets_begin(), p.facets_end(), p);
930
931 edgeIntersectionsCGAL
932 (
933 tree,
934 surf1,
935 surf1.points(),
936 edgeCuts1
937 );
938 }
939 // See which edges of 2 intersect 1
940 {
941 Info<< "Intersect surface 2 edges with surface 1:" << nl;
942 Info<< " constructing CGAL surface ..." << endl;
944 PolyhedronReader(surf1, p);
945
946 Info<< " constructing CGAL tree ..." << endl;
947 const Tree tree(p.facets_begin(), p.facets_end(), p);
948
949 edgeIntersectionsCGAL
950 (
951 tree,
952 surf2,
953 surf2.points(),
954 edgeCuts2
955 );
956 }
957 Info<< endl;
958 }
959 else
960 {
961 const scalarField surf1PointTol
962 (
963 max(1e-8*edgeIntersections::minEdgeLength(surf1), SMALL)
964 );
965 const scalarField surf2PointTol
966 (
967 max(1e-8*edgeIntersections::minEdgeLength(surf2), SMALL)
968 );
969
970
971 Random rndGen(0);
972
973 labelPair cuts1;
974 labelPair cuts2;
975
976 for (label iter = 0; iter < 10; iter++)
977 {
978 // Find intersections of surf1 edges with surf2 triangles
979 cuts1 = edgeIntersectionsAndShuffleCGAL
980 (
981 rndGen,
982 surf2,
983 surf1PointTol,
984 surf1,
985 edgeCuts1
986 );
987
988 // Find intersections of surf2 edges with surf1 triangles
989 cuts2 = edgeIntersectionsAndShuffleCGAL
990 (
991 rndGen,
992 surf1,
993 surf2PointTol,
994 surf2,
995 edgeCuts2
996 );
997
998 if (cuts1.second() + cuts2.second() == 0)
999 {
1000 break;
1001 }
1002 }
1003 }
1004}
1005
1006
1007void calcEdgeCutsBitsCGAL
1008(
1009 triSurface& surf1,
1010 triSurface& surf2,
1011 const bool perturb,
1012 edgeIntersections& edgeCuts1,
1013 edgeIntersections& edgeCuts2
1014)
1015{
1016 label nZones1 = 1;
1017 labelList zone1;
1018 {
1019 labelHashSet orientationEdge(surf1.size()/1000);
1020 PatchTools::checkOrientation(surf1, false, &orientationEdge);
1021 nZones1 = PatchTools::markZones(surf1, orientationEdge, zone1);
1022
1023 Info<< "Surface triangles " << surf1.size()
1024 << " number of zones (orientation or non-manifold):"
1025 << nZones1 << endl;
1026 }
1027
1028 label nZones2 = 1;
1029 labelList zone2;
1030 {
1031 labelHashSet orientationEdge(surf2.size()/1000);
1032 PatchTools::checkOrientation(surf2, false, &orientationEdge);
1033 nZones2 = PatchTools::markZones(surf2, orientationEdge, zone2);
1034
1035 Info<< "Surface triangles " << surf2.size()
1036 << " number of zones (orientation or non-manifold):"
1037 << nZones2 << endl;
1038 }
1039
1040
1041 if (nZones1 == 1 && nZones2 == 1)
1042 {
1043 calcEdgeCutsCGAL
1044 (
1045 surf1,
1046 surf2,
1047 perturb,
1048 edgeCuts1,
1049 edgeCuts2
1050 );
1051 }
1052 else
1053 {
1054 edgeCuts1 = edgeIntersections
1055 (
1056 List<List<pointIndexHit>>(surf1.nEdges()),
1057 labelListList(surf1.nEdges())
1058 );
1059 edgeCuts2 = edgeIntersections
1060 (
1061 List<List<pointIndexHit>>(surf2.nEdges()),
1062 labelListList(surf2.nEdges())
1063 );
1064
1065
1066 for (label zone1I = 0; zone1I < nZones1; zone1I++)
1067 {
1068 // Generate sub surface for zone1I
1069 boolList includeMap1(surf1.size(), false);
1070
1071 forAll(zone1, faceI)
1072 {
1073 if (zone1[faceI] == zone1I)
1074 {
1075 includeMap1[faceI] = true;
1076 }
1077 }
1078
1079 // Subset. Map from local points on subset to local points on
1080 // original
1081 labelList pointMap1;
1082 labelList faceMap1;
1083 triSurface subSurf1
1084 (
1085 surf1.subsetMesh
1086 (
1087 includeMap1,
1088 pointMap1,
1089 faceMap1
1090 )
1091 );
1092
1093 // Split non-manifold points; update pointMap
1094 dupNonManifoldPoints(subSurf1, pointMap1);
1095
1096 const boundBox subBb1
1097 (
1098 subSurf1.points(),
1099 subSurf1.meshPoints(),
1100 false
1101 );
1102
1103 const labelList edgeMap1
1104 (
1105 matchEdges
1106 (
1107 surf1,
1108 subSurf1,
1109 pointMap1
1110 )
1111 );
1112
1113
1114 for (label zone2I = 0; zone2I < nZones2; zone2I++)
1115 {
1116 // Generate sub surface for zone2I
1117 boolList includeMap2(surf2.size(), false);
1118
1119 forAll(zone2, faceI)
1120 {
1121 if (zone2[faceI] == zone2I)
1122 {
1123 includeMap2[faceI] = true;
1124 }
1125 }
1126
1127 labelList pointMap2;
1128 labelList faceMap2;
1129 triSurface subSurf2
1130 (
1131 surf2.subsetMesh
1132 (
1133 includeMap2,
1134 pointMap2,
1135 faceMap2
1136 )
1137 );
1138
1139
1140 const boundBox subBb2
1141 (
1142 subSurf2.points(),
1143 subSurf2.meshPoints(),
1144 false
1145 );
1146
1147 // Short-circuit expensive calculations
1148 if (!subBb2.overlaps(subBb1))
1149 {
1150 continue;
1151 }
1152
1153
1154 // Split non-manifold points; update pointMap
1155 dupNonManifoldPoints(subSurf2, pointMap2);
1156
1157 const labelList edgeMap2
1158 (
1159 matchEdges
1160 (
1161 surf2,
1162 subSurf2,
1163 pointMap2
1164 )
1165 );
1166
1167
1168 // Do cuts
1169 edgeIntersections subEdgeCuts1;
1170 edgeIntersections subEdgeCuts2;
1171 calcEdgeCutsCGAL
1172 (
1173 subSurf1,
1174 subSurf2,
1175 perturb,
1176 subEdgeCuts1,
1177 subEdgeCuts2
1178 );
1179
1180 // Move original surface
1181 {
1182 pointField points2(surf2.points());
1183 forAll(pointMap2, i)
1184 {
1185 label subMeshPointI = subSurf2.meshPoints()[i];
1186 label meshPointI = surf2.meshPoints()[pointMap2[i]];
1187 points2[meshPointI] = subSurf2.points()[subMeshPointI];
1188 }
1189 surf2.movePoints(points2);
1190 }
1191
1192 // Insert into main structure
1193 edgeCuts1.merge(subEdgeCuts1, edgeMap1, faceMap2);
1194 edgeCuts2.merge(subEdgeCuts2, edgeMap2, faceMap1);
1195 }
1196
1197
1198 // Move original surface
1199 {
1200 pointField points1(surf1.points());
1201 forAll(pointMap1, i)
1202 {
1203 label subMeshPointI = subSurf1.meshPoints()[i];
1204 label meshPointI = surf1.meshPoints()[pointMap1[i]];
1205 points1[meshPointI] = subSurf1.points()[subMeshPointI];
1206 }
1207 surf1.movePoints(points1);
1208 }
1209 }
1210 }
1211}
1212
1213#endif // NO_CGAL
1214
1215
1216//void calcFeaturePoints(const pointField& points, const edgeList& edges)
1217//{
1218// edgeMesh eMesh(points, edges);
1219//
1220// const labelListList& pointEdges = eMesh.pointEdges();
1221//
1222//
1223// // Get total number of feature points
1224// label nFeaturePoints = 0;
1225// forAll(pointEdges, pI)
1226// {
1227// const labelList& pEdges = pointEdges[pI];
1228//
1229// if (pEdges.size() == 1)
1230// {
1231// nFeaturePoints++;
1232// }
1233// }
1234//
1235//
1236// // Calculate addressing from feature point to cut point and cut edge
1237// labelList featurePointToCutPoint(nFeaturePoints);
1238// labelList featurePointToCutEdge(nFeaturePoints);
1239//
1240// label nFeatPts = 0;
1241// forAll(pointEdges, pI)
1242// {
1243// const labelList& pEdges = pointEdges[pI];
1244//
1245// if (pEdges.size() == 1)
1246// {
1247// featurePointToCutPoint[nFeatPts] = pI;
1248// featurePointToCutEdge[nFeatPts] = pEdges[0];
1249// nFeatPts++;
1250// }
1251// }
1252//
1253//
1254//
1255// label concaveStart = 0;
1256// label mixedStart = 0;
1257// label nonFeatureStart = nFeaturePoints;
1258//
1259//
1260// labelListList featurePointNormals(nFeaturePoints);
1261// labelListList featurePointEdges(nFeaturePoints);
1262// labelList regionEdges;
1263//}
1264
1265
1267(
1268 const IOobject& io,
1269 const booleanSurface::booleanOpType action,
1270 const bool surf1Baffle,
1271 const bool surf2Baffle,
1272 const bool invertedSpace,
1273 const triSurface& surf1,
1274 const edgeIntersections& edgeCuts1,
1275 const triSurface& surf2,
1276 const edgeIntersections& edgeCuts2
1277)
1278{
1279 // Determine intersection edges from the edge cuts
1281 (
1282 surf1,
1283 edgeCuts1,
1284 surf2,
1285 edgeCuts2
1286 );
1287
1288 label nFeatEds = inter.cutEdges().size();
1289
1290 DynamicList<vector> normals(2*nFeatEds);
1291 vectorField edgeDirections(nFeatEds, Zero);
1293 (
1294 2*nFeatEds
1295 );
1296 List<DynamicList<label>> edgeNormals(nFeatEds);
1297 List<DynamicList<label>> normalDirections(nFeatEds);
1298
1299
1300 const triSurface& s1 = surf1;
1301 const triSurface& s2 = surf2;
1302
1303 forAllConstIters(inter.facePairToEdgeId(), iter)
1304 {
1305 const labelPair& facePair = iter.key();
1306 const label cutEdgeI = iter.val();
1307
1308 const edge& fE = inter.cutEdges()[cutEdgeI];
1309
1310 const vector& norm1 = s1.faceNormals()[facePair.first()];
1311 const vector& norm2 = s2.faceNormals()[facePair.second()];
1312
1313 DynamicList<label>& eNormals = edgeNormals[cutEdgeI];
1314 DynamicList<label>& nDirections = normalDirections[cutEdgeI];
1315
1316 edgeDirections[cutEdgeI] = fE.vec(inter.cutPoints());
1317
1318 normals.append(norm1);
1319 eNormals.append(normals.size() - 1);
1320
1321 if (surf1Baffle)
1322 {
1323 normalVolumeTypes.append(extendedFeatureEdgeMesh::BOTH);
1324
1325 nDirections.append(1);
1326 }
1327 else
1328 {
1329 normalVolumeTypes.append(extendedFeatureEdgeMesh::INSIDE);
1330 nDirections.append
1331 (
1332 calcNormalDirection
1333 (
1334 norm1,
1335 norm2,
1336 edgeDirections[cutEdgeI],
1337 s1[facePair.first()].centre(s1.points()),
1338 inter.cutPoints()[fE.start()]
1339 )
1340 );
1341 }
1342
1343 normals.append(norm2);
1344 eNormals.append(normals.size() - 1);
1345
1346 if (surf2Baffle)
1347 {
1348 normalVolumeTypes.append(extendedFeatureEdgeMesh::BOTH);
1349
1350 nDirections.append(1);
1351 }
1352 else
1353 {
1354 normalVolumeTypes.append(extendedFeatureEdgeMesh::INSIDE);
1355
1356 nDirections.append
1357 (
1358 calcNormalDirection
1359 (
1360 norm2,
1361 norm1,
1362 edgeDirections[cutEdgeI],
1363 s2[facePair.second()].centre(s2.points()),
1364 inter.cutPoints()[fE.start()]
1365 )
1366 );
1367 }
1368
1369
1370 if (surf1Baffle)
1371 {
1372 normals.append(norm2);
1373
1374 if (surf2Baffle)
1375 {
1376 normalVolumeTypes.append(extendedFeatureEdgeMesh::BOTH);
1377
1378 nDirections.append(1);
1379 }
1380 else
1381 {
1382 normalVolumeTypes.append(extendedFeatureEdgeMesh::INSIDE);
1383
1384 nDirections.append
1385 (
1386 calcNormalDirection
1387 (
1388 norm2,
1389 norm1,
1390 edgeDirections[cutEdgeI],
1391 s2[facePair.second()].centre(s2.points()),
1392 inter.cutPoints()[fE.start()]
1393 )
1394 );
1395 }
1396
1397 eNormals.append(normals.size() - 1);
1398 }
1399
1400 if (surf2Baffle)
1401 {
1402 normals.append(norm1);
1403
1404 if (surf1Baffle)
1405 {
1406 normalVolumeTypes.append(extendedFeatureEdgeMesh::BOTH);
1407
1408 nDirections.append(1);
1409 }
1410 else
1411 {
1412 normalVolumeTypes.append(extendedFeatureEdgeMesh::INSIDE);
1413
1414 nDirections.append
1415 (
1416 calcNormalDirection
1417 (
1418 norm1,
1419 norm2,
1420 edgeDirections[cutEdgeI],
1421 s1[facePair.first()].centre(s1.points()),
1422 inter.cutPoints()[fE.start()]
1423 )
1424 );
1425 }
1426
1427 eNormals.append(normals.size() - 1);
1428 }
1429 }
1430
1431
1432 label internalStart = -1;
1433 label nIntOrExt = 0;
1434 label nFlat = 0;
1435 label nOpen = 0;
1436 // label nMultiple = 0;
1437
1438 forAll(edgeNormals, eI)
1439 {
1440 label nEdNorms = edgeNormals[eI].size();
1441
1442 if (nEdNorms == 1)
1443 {
1444 nOpen++;
1445 }
1446 else if (nEdNorms == 2)
1447 {
1448 const vector& n0(normals[edgeNormals[eI][0]]);
1449 const vector& n1(normals[edgeNormals[eI][1]]);
1450
1452 {
1453 nFlat++;
1454 }
1455 else
1456 {
1457 nIntOrExt++;
1458 }
1459 }
1460 // else if (nEdNorms > 2)
1461 // {
1462 // ++nMultiple;
1463 // }
1464 }
1465
1466 if (action == booleanSurface::UNION)
1467 {
1468 if (!invertedSpace)
1469 {
1470 // All edges are internal
1471 internalStart = 0;
1472 }
1473 else
1474 {
1475 // All edges are external
1476 internalStart = nIntOrExt;
1477 }
1478 }
1479 else if (action == booleanSurface::INTERSECTION)
1480 {
1481 if (!invertedSpace)
1482 {
1483 // All edges are external
1484 internalStart = nIntOrExt;
1485 }
1486 else
1487 {
1488 // All edges are internal
1489 internalStart = 0;
1490 }
1491 }
1492 else if (action == booleanSurface::DIFFERENCE)
1493 {
1494 // All edges are external
1495 internalStart = nIntOrExt;
1496 }
1497 else
1498 {
1500 << "Unsupported booleanSurface:booleanOpType and space "
1501 << action << " " << invertedSpace
1502 << abort(FatalError);
1503 }
1504
1505 // There are no feature points supported by surfaceIntersection
1506 // Flat, open or multiple edges are assumed to be impossible
1507 // Region edges are not explicitly supported by surfaceIntersection
1508
1509 vectorField normalsTmp(normals);
1510 List<extendedFeatureEdgeMesh::sideVolumeType> normalVolumeTypesTmp
1511 (
1512 normalVolumeTypes
1513 );
1514 labelListList edgeNormalsTmp(edgeNormals.size());
1515 forAll(edgeNormalsTmp, i)
1516 {
1517 edgeNormalsTmp[i] = edgeNormals[i];
1518 }
1519 labelListList normalDirectionsTmp(normalDirections.size());
1520 forAll(normalDirectionsTmp, i)
1521 {
1522 normalDirectionsTmp[i] = normalDirections[i];
1523 }
1524
1525 //calcFeaturePoints(inter.cutPoints(), inter.cutEdges());
1526
1528 (
1529 io,
1530 inter.cutPoints(),
1531 inter.cutEdges(),
1532
1533 0, // concaveStart,
1534 0, // mixedStart,
1535 0, // nonFeatureStart,
1536
1537 internalStart, // internalStart,
1538 nIntOrExt, // flatStart,
1539 nIntOrExt + nFlat, // openStart,
1540 nIntOrExt + nFlat + nOpen, // multipleStart,
1541
1542 normalsTmp,
1543 normalVolumeTypesTmp,
1544 edgeDirections,
1545 normalDirectionsTmp,
1546 edgeNormalsTmp,
1547
1548 labelListList(), // featurePointNormals,
1549 labelListList(), // featurePointEdges,
1550 labelList() // regionEdges
1551 );
1552}
1553
1554int main(int argc, char *argv[])
1555{
1557 (
1558 "Generates a extendedFeatureEdgeMesh for the interface created by"
1559 " a boolean operation on two surfaces."
1560 #ifndef NO_CGAL
1561 " [Compiled with CGAL]"
1562 #else
1563 " [Compiled without CGAL]"
1564 #endif
1565 );
1566
1569 (
1570 "action",
1571 "One of (intersection | union | difference)"
1572 );
1573 argList::addArgument("surface1", "The input surface file 1");
1574 argList::addArgument("surface2", "The input surface file 2");
1575
1577 (
1578 "scale",
1579 "factor",
1580 "Geometry scaling factor (both surfaces)"
1581 );
1583 (
1584 "surf1Baffle",
1585 "Mark surface 1 as a baffle"
1586 );
1587
1589 (
1590 "surf2Baffle",
1591 "Mark surface 2 as a baffle"
1592 );
1593
1595 (
1596 "perturb",
1597 "Perturb surface points to escape degenerate intersections"
1598 );
1599
1601 (
1602 "no-cgal",
1603 #ifndef NO_CGAL
1604 "Do not use CGAL algorithms"
1605 #else
1606 "Ignored, compiled without CGAL"
1607 #endif
1608 );
1609
1611 (
1612 "invertedSpace",
1613 "Do the surfaces have inverted space orientation, "
1614 "i.e. a point at infinity is considered inside. "
1615 "This is only sensible for union and intersection."
1616 );
1617
1619 (
1620 "trim",
1621 "((surface1 volumeType) .. (surfaceN volumeType))",
1622 "Trim resulting intersection with additional surfaces;"
1623 " volumeType is 'inside' (keep (parts of) edges that are inside)"
1624 ", 'outside' (keep (parts of) edges that are outside) or"
1625 " 'mixed' (keep all)"
1626 );
1627
1628 #include "setRootCase.H"
1629 #include "createTime.H"
1630
1631 const word action(args[1]);
1632
1633 const Enum<booleanSurface::booleanOpType> validActions
1634 {
1635 { booleanSurface::INTERSECTION, "intersection" },
1636 { booleanSurface::UNION, "union" },
1637 { booleanSurface::DIFFERENCE, "difference" }
1638 };
1639
1640 if (!validActions.found(action))
1641 {
1643 << "Unsupported action " << action << endl
1644 << "Supported actions:" << validActions << nl
1645 << abort(FatalError);
1646 }
1647
1648
1649 List<Pair<word>> surfaceAndSide;
1650 if (args.readIfPresent("trim", surfaceAndSide))
1651 {
1652 Info<< "Trimming edges with " << surfaceAndSide << endl;
1653 }
1654
1655
1656 // Scale factor for both surfaces:
1657 const scalar scaleFactor = args.getOrDefault<scalar>("scale", -1);
1658
1659 const word surf1Name(args[2]);
1660 Info<< "Reading surface " << surf1Name << endl;
1661 triSurfaceMesh surf1
1662 (
1663 IOobject
1664 (
1665 surf1Name,
1666 runTime.constant(),
1668 runTime
1669 )
1670 );
1671 if (scaleFactor > 0)
1672 {
1673 Info<< "Scaling : " << scaleFactor << nl;
1674 surf1.scalePoints(scaleFactor);
1675 }
1676
1677 Info<< surf1Name << " statistics:" << endl;
1678 surf1.writeStats(Info);
1679 Info<< endl;
1680
1681 const word surf2Name(args[3]);
1682 Info<< "Reading surface " << surf2Name << endl;
1683 triSurfaceMesh surf2
1684 (
1685 IOobject
1686 (
1687 surf2Name,
1688 runTime.constant(),
1690 runTime
1691 )
1692 );
1693 if (scaleFactor > 0)
1694 {
1695 Info<< "Scaling : " << scaleFactor << nl;
1696 surf2.scalePoints(scaleFactor);
1697 }
1698
1699 Info<< surf2Name << " statistics:" << endl;
1700 surf2.writeStats(Info);
1701 Info<< endl;
1702
1703 const bool surf1Baffle = args.found("surf1Baffle");
1704 const bool surf2Baffle = args.found("surf2Baffle");
1705
1706 edgeIntersections edgeCuts1;
1707 edgeIntersections edgeCuts2;
1708
1709 const bool invertedSpace = args.found("invertedSpace");
1710
1711 if (invertedSpace && validActions[action] == booleanSurface::DIFFERENCE)
1712 {
1714 << "Inverted space only makes sense for union or intersection."
1715 << exit(FatalError);
1716 }
1717
1718
1719 // Calculate where edges are cut by the other surface
1720#ifndef NO_CGAL
1721 if (!args.found("no-cgal"))
1722 {
1723 calcEdgeCutsBitsCGAL
1724 (
1725 surf1,
1726 surf2,
1727 args.found("perturb"),
1728 edgeCuts1,
1729 edgeCuts2
1730 );
1731 }
1732 else
1733#endif // NO_CGAL
1734 {
1735 calcEdgeCuts
1736 (
1737 surf1,
1738 surf2,
1739 args.found("perturb"),
1740 edgeCuts1,
1741 edgeCuts2
1742 );
1743 }
1744
1745 const fileName sFeatFileName
1746 (
1747 fileName::stem(surf1Name)
1748 + "_"
1749 + fileName::stem(surf2Name)
1750 + "_"
1751 + action
1752 );
1753
1755 (
1756 createEdgeMesh
1757 (
1758 IOobject
1759 (
1760 sFeatFileName + ".extendedFeatureEdgeMesh",
1761 runTime.constant(),
1762 "extendedFeatureEdgeMesh",
1763 runTime,
1766 ),
1768 surf1Baffle,
1769 surf2Baffle,
1770 invertedSpace,
1771 surf1,
1772 edgeCuts1,
1773 surf2,
1774 edgeCuts2
1775 )
1776 );
1777
1778
1779 // Trim intersections
1780 forAll(surfaceAndSide, trimI)
1781 {
1782 const word& trimName = surfaceAndSide[trimI].first();
1783 const volumeType trimType
1784 (
1785 volumeType::names[surfaceAndSide[trimI].second()]
1786 );
1787
1788 Info<< "Reading trim surface " << trimName << endl;
1789 const triSurfaceMesh trimSurf
1790 (
1791 IOobject
1792 (
1793 trimName,
1794 runTime.constant(),
1796 runTime,
1799 )
1800 );
1801
1802 Info<< trimName << " statistics:" << endl;
1803 trimSurf.writeStats(Info);
1804 Info<< endl;
1805
1806 labelList pointMap;
1807 labelList edgeMap;
1808 feMeshPtr().trim
1809 (
1810 trimSurf,
1811 trimType,
1812 pointMap,
1813 edgeMap
1814 );
1815 }
1816
1817
1818 const extendedFeatureEdgeMesh& feMesh = feMeshPtr();
1819
1820 feMesh.writeStats(Info);
1821 Info<< nl << "Writing extendedFeatureEdgeMesh to "
1822 << feMesh.objectRelPath() << nl
1823 << endl;
1824 feMesh.write();
1825 feMesh.writeObj(feMesh.path()/sFeatFileName);
1826
1827 {
1828 // Write a featureEdgeMesh for backwards compatibility
1829 featureEdgeMesh bfeMesh
1830 (
1831 IOobject
1832 (
1833 sFeatFileName + ".eMesh", // name
1834 runTime.constant(), // instance
1836 runTime, // registry
1840 ),
1841 feMesh.points(),
1842 feMesh.edges()
1843 );
1844
1845 Info<< nl << "Writing featureEdgeMesh to "
1846 << bfeMesh.objectRelPath() << endl;
1847
1848 bfeMesh.regIOobject::write();
1849 }
1850
1851 Info << nl << "End\n" << endl;
1852
1853 return 0;
1854}
1855
1856
1857// ************************************************************************* //
CGAL data structures used for triSurface handling.
CGAL::Segment_3< K > Segment
CGAL::Polyhedron_3< K, My_items > Polyhedron
CGAL::Point_3< K > Point
Dynamically sized Field. Similar to DynamicList, but inheriting from a Field instead of a List.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
bool found(const word &key) const
Same as contains().
Definition Enum.H:480
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
fileName objectRelPath() const
The object path relative to the case.
Definition IOobject.C:581
fileName path() const
The complete path for the object (with instance, local,...).
Definition IOobject.C:500
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
const T & first() const noexcept
Access the first element.
Definition Pair.H:137
const T & second() const noexcept
Access the second element.
Definition Pair.H:147
static bool checkOrientation(const PrimitivePatch< FaceList, PointField > &, const bool report=false, labelHashSet *marked=0)
Check for orientation issues.
static label markZones(const PrimitivePatch< FaceList, PointField > &, const BoolListType &borderEdge, labelList &faceZone)
Size and fills faceZone with zone of face.
label nEdges() const
Number of edges in patch.
label nPoints() const
Number of points supporting patch faces.
const labelListList & pointEdges() const
Return point-edge addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
Random number generator.
Definition Random.H:56
T & first()
Access first element of the list, position [0].
Definition UList.H:957
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition argList.C:366
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition argList.C:389
static void noParallel()
Remove the parallel options.
Definition argList.C:599
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition argList.C:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
static const Enum< booleanOpType > booleanOpTypeNames
booleanOpType
Enumeration listing the possible volume operator types.
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
Holder of intersections of edges of a surface with another surface. Optionally shuffles around points...
const labelListList & classification() const
For every intersection the classification status.
void merge(const edgeIntersections &, const labelList &edgeMap, const labelList &faceMap, const bool merge=true)
Merge (or override) edge intersection for a subset.
static scalarField minEdgeLength(const triSurface &surf)
Calculate min edge length for every surface point.
label removeDegenerates(const label nIters, const triSurface &surf1, const triSurfaceSearch &query2, const scalarField &surf1PointTol, pointField &points1)
Resolve ties. Shuffles points so all edge - face intersections.
const pointField & points() const noexcept
Return points.
Definition edgeMeshI.H:92
const edgeList & edges() const noexcept
Return edges.
Definition edgeMeshI.H:98
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
vector vec(const UList< point > &pts) const
Return the vector (from first to second).
Definition edgeI.H:403
linePointRef line(const UList< point > &pts) const
Return edge line.
Definition edgeI.H:463
label start() const noexcept
The start (first) vertex label.
Definition edge.H:155
void writeObj(const fileName &prefix) const
Write all components of the extendedEdgeMesh as obj files.
virtual void writeStats(Ostream &os) const
Dump some information.
static scalar cosNormalAngleTol_
Angular closeness tolerance for treating normals as the same.
virtual bool write(const bool writeOnProc=true) const
Give precedence to the regIOobject write.
A class for handling file names.
Definition fileName.H:75
static std::string stem(const std::string &str)
Return the basename, without extension.
Definition fileName.C:391
A triFace with additional (region) index.
Definition labelledTri.H:56
Basic surface-surface intersection description. Constructed from two surfaces it creates a descriptio...
IOoject and searching on triSurface.
static word meshSubDir
Return the mesh sub-directory name (usually "triSurface").
Helper class to search on triSurface.
Triangulated surface description with patch information.
Definition triSurface.H:74
triSurface subsetMesh(const UList< bool > &include, labelList &pointMap, labelList &faceMap) const
Return a new surface subsetted on the selected faces.
Definition triSurface.C:873
void write(Ostream &os) const
Write to Ostream in simple OpenFOAM format.
void writeStats(Ostream &os) const
Write some statistics.
virtual void movePoints(const pointField &pts)
Move points.
Definition triSurface.C:606
virtual void scalePoints(const scalar scaleFactor)
Scale points. A non-positive factor is ignored.
Definition triSurface.C:632
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition volumeType.H:56
static const Enum< volumeType::type > names
Names for the classification enumeration.
Definition volumeType.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
const pointField & points
label nPoints
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const char * end
Definition SVGTools.H:223
const dimensionedScalar mp
Proton mass.
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
UIndirectList< label > labelUIndList
UIndirectList of labels.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void cross(FieldField< Field1, typename crossProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};constexpr label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< DynamicList< face > > pFaces[nBCs]
labelList f(nPoints)
Tree tree(triangles.begin(), triangles.end())
Foam::argList args(argc, argv)
volScalarField & b
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Random rndGen