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