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