booleanSurface.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 | Copyright (C) 2015 OpenCFD Ltd.
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 \*---------------------------------------------------------------------------*/
25 
26 #include "booleanSurface.H"
27 #include "intersectedSurface.H"
28 #include "orientedSurface.H"
29 #include "triSurfaceSearch.H"
30 #include "OFstream.H"
31 #include "treeBoundBox.H"
32 #include "meshTools.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 defineTypeNameAndDebug(booleanSurface, 0);
39 
40 template<>
41 const char* Foam::NamedEnum
42 <
44  4
45 >::names[] =
46 {
47  "union",
48  "intersection",
49  "difference",
50  "all"
51 };
52 }
53 
56 
57 
58 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59 
60 // Check whether at least one of faces connected to the intersection has been
61 // marked.
63 (
64  const intersectedSurface& surf,
65  const labelList& faceZone,
66  const label includedFace
67 )
68 {
69  forAll(surf.intersectionEdges(), intEdgeI)
70  {
71  label edgeI = surf.intersectionEdges()[intEdgeI];
72 
73  const labelList& myFaces = surf.edgeFaces()[edgeI];
74 
75  bool usesIncluded = false;
76 
77  forAll(myFaces, myFaceI)
78  {
79  if (faceZone[myFaces[myFaceI]] == faceZone[includedFace])
80  {
81  usesIncluded = true;
82 
83  break;
84  }
85  }
86 
87  if (!usesIncluded)
88  {
90  << "None of the faces reachable from face " << includedFace
91  << " connects to the intersection."
92  << exit(FatalError);
93  }
94  }
95 }
96 
97 
98 // Linear lookup
100 (
101  const labelList& elems,
102  const label elem
103 )
104 {
105  forAll(elems, elemI)
106  {
107  if (elems[elemI] == elem)
108  {
109  return elemI;
110  }
111  }
112  return -1;
113 }
114 
115 
117 (
118  const edgeList& edges,
119  const labelList& edgeLabels,
120  const edge& e
121 )
122 {
123  forAll(edgeLabels, edgeLabelI)
124  {
125  if (edges[edgeLabels[edgeLabelI]] == e)
126  {
127  return edgeLabels[edgeLabelI];
128  }
129  }
131  << "Cannot find edge " << e << " in edges " << edgeLabels
132  << abort(FatalError);
133 
134  return -1;
135 }
136 
137 
138 // Generate combined patchList (returned). Sets patchMap to map from surf2
139 // region numbers into combined patchList
141 (
142  const triSurface& surf1,
143  const triSurface& surf2,
144  labelList& patchMap2
145 )
146 {
147  // Size too big.
148  geometricSurfacePatchList combinedPatches
149  (
150  surf1.patches().size()
151  + surf2.patches().size()
152  );
153 
154  // Copy all patches of surf1
155  label combinedPatchI = 0;
156  forAll(surf1.patches(), patchI)
157  {
158  combinedPatches[combinedPatchI++] = surf1.patches()[patchI];
159  }
160 
161  // (inefficiently) add unique patches from surf2
162  patchMap2.setSize(surf2.patches().size());
163 
164  forAll(surf2.patches(), patch2I)
165  {
166  label index = -1;
167 
168  forAll(surf1.patches(), patch1I)
169  {
170  if (surf1.patches()[patch1I] == surf2.patches()[patch2I])
171  {
172  index = patch1I;
173 
174  break;
175  }
176  }
177 
178  if (index == -1)
179  {
180  combinedPatches[combinedPatchI] = surf2.patches()[patch2I];
181  patchMap2[patch2I] = combinedPatchI;
182  combinedPatchI++;
183  }
184  else
185  {
186  patchMap2[patch2I] = index;
187  }
188  }
189 
190  combinedPatches.setSize(combinedPatchI);
191 
192  return combinedPatches;
193 }
194 
195 
197 (
198  const triSurface& surf,
199  const label prevVert0,
200  const label prevFaceI,
201  const label prevState,
202  const label edgeI,
203  labelList& side
204 )
205 {
206  const labelList& eFaces = surf.sortedEdgeFaces()[edgeI];
207 
208  // Simple case. Propagate side.
209  if (eFaces.size() == 2)
210  {
211  forAll(eFaces, eFaceI)
212  {
213  propagateSide
214  (
215  surf,
216  prevState,
217  eFaces[eFaceI],
218  side
219  );
220  }
221  }
222 
223 
224  if (((eFaces.size() % 2) == 1) && (eFaces.size() != 1))
225  {
227  << "Don't know how to handle edges with odd number of faces"
228  << endl
229  << "edge:" << edgeI << " vertices:" << surf.edges()[edgeI]
230  << " coming from face:" << prevFaceI
231  << " edgeFaces:" << eFaces << abort(FatalError);
232  }
233 
234 
235  // Get position of face in edgeFaces
236  label ind = index(eFaces, prevFaceI);
237 
238  // Determine orientation of faces around edge prevVert0
239  // (might be opposite of edge)
240  const edge& e = surf.edges()[edgeI];
241 
242  // Get next face to include
243  label nextInd;
244  label prevInd;
245 
246  if (e.start() == prevVert0)
247  {
248  // Edge (and hence eFaces) in same order as prevVert0.
249  // Take next face from sorted list
250  nextInd = eFaces.fcIndex(ind);
251  prevInd = eFaces.rcIndex(ind);
252  }
253  else
254  {
255  // Take previous face from sorted neighbours
256  nextInd = eFaces.rcIndex(ind);
257  prevInd = eFaces.fcIndex(ind);
258  }
259 
260 
261  if (prevState == OUTSIDE)
262  {
263  // Coming from outside. nextInd is outside, rest is inside.
264 
265  forAll(eFaces, eFaceI)
266  {
267  if (eFaceI != ind)
268  {
269  label nextState;
270 
271  if (eFaceI == nextInd)
272  {
273  nextState = OUTSIDE;
274  }
275  else
276  {
277  nextState = INSIDE;
278  }
279 
280  propagateSide
281  (
282  surf,
283  nextState,
284  eFaces[eFaceI],
285  side
286  );
287  }
288  }
289  }
290  else
291  {
292  // Coming from inside. prevInd is inside as well, rest is outside.
293 
294  forAll(eFaces, eFaceI)
295  {
296  if (eFaceI != ind)
297  {
298  label nextState;
299 
300  if (eFaceI == prevInd)
301  {
302  nextState = INSIDE;
303  }
304  else
305  {
306  nextState = OUTSIDE;
307  }
308 
309  propagateSide
310  (
311  surf,
312  nextState,
313  eFaces[eFaceI],
314  side
315  );
316  }
317  }
318  }
319 }
320 
321 
322 // Face-edge walk. Determines inside/outside for all faces connected to an edge.
324 (
325  const triSurface& surf,
326  const label prevState,
327  const label faceI,
328  labelList& side
329 )
330 {
331  if (side[faceI] == UNVISITED)
332  {
333  side[faceI] = prevState;
334 
335  const labelledTri& tri = surf.localFaces()[faceI];
336 
337  // Get copy of face labels
338  label a = tri[0];
339  label b = tri[1];
340  label c = tri[2];
341 
342  // Go and visit my edges' face-neighbours.
343 
344  const labelList& myEdges = surf.faceEdges()[faceI];
345 
346  label edgeAB = findEdge(surf.edges(), myEdges, edge(a, b));
347 
348  propagateEdgeSide
349  (
350  surf,
351  a,
352  faceI,
353  prevState,
354  edgeAB,
355  side
356  );
357 
358  label edgeBC = findEdge(surf.edges(), myEdges, edge(b, c));
359 
360  propagateEdgeSide
361  (
362  surf,
363  b,
364  faceI,
365  prevState,
366  edgeBC,
367  side
368  );
369 
370  label edgeCA = findEdge(surf.edges(), myEdges, edge(c, a));
371 
372  propagateEdgeSide
373  (
374  surf,
375  c,
376  faceI,
377  prevState,
378  edgeCA,
379  side
380  );
381  }
382 }
383 
384 
385 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
386 
387 // Null constructor
389 :
390  triSurface()
391 {}
392 
393 
394 // Construct from surfaces and face to include for every surface
396 (
397  const triSurface& surf1,
398  const triSurface& surf2,
399  const surfaceIntersection& inter,
400  const label includeFace1,
401  const label includeFace2
402 )
403 :
404  triSurface(),
405  faceMap_()
406 {
407  if (debug)
408  {
409  Pout<< "booleanSurface : Generating intersected surface for surf1"
410  << endl;
411  }
412 
413  // Add intersection to surface1 (retriangulates cut faces)
414  intersectedSurface cutSurf1(surf1, true, inter);
415 
416 
417  if (debug)
418  {
419  Pout<< "booleanSurface : Generated cutSurf1: " << endl;
420  cutSurf1.writeStats(Pout);
421 
422  Pout<< "Writing to file cutSurf1.obj" << endl;
423  cutSurf1.write("cutSurf1.obj");
424  }
425 
426  if (debug)
427  {
428  Pout<< "booleanSurface : Generating intersected surface for surf2"
429  << endl;
430  }
431 
432  // Add intersection to surface2
433  intersectedSurface cutSurf2(surf2, false, inter);
434 
435  if (debug)
436  {
437  Pout<< "booleanSurface : Generated cutSurf2: " << endl;
438  cutSurf2.writeStats(Pout);
439 
440  Pout<< "Writing to file cutSurf2.obj" << endl;
441  cutSurf2.write("cutSurf2.obj");
442  }
443 
444 
445  // Find (first) face of cutSurf1 that originates from includeFace1
446  label cutSurf1FaceI = index(cutSurf1.faceMap(), includeFace1);
447 
448  if (debug)
449  {
450  Pout<< "cutSurf1 : starting to fill from face:" << cutSurf1FaceI
451  << endl;
452  }
453 
454  if (cutSurf1FaceI == -1)
455  {
457  << "Did not find face with label " << includeFace1
458  << " in intersectedSurface."
459  << exit(FatalError);
460  }
461 
462  // Find (first) face of cutSurf2 that originates from includeFace1
463  label cutSurf2FaceI = index(cutSurf2.faceMap(), includeFace2);
464 
465  if (debug)
466  {
467  Pout<< "cutSurf2 : starting to fill from face:" << cutSurf2FaceI
468  << endl;
469  }
470  if (cutSurf2FaceI == -1)
471  {
473  << "Did not find face with label " << includeFace2
474  << " in intersectedSurface."
475  << exit(FatalError);
476  }
477 
478 
479  //
480  // Mark faces of cutSurf1 that need to be kept by walking from includeFace1
481  // without crossing any edges of the intersection.
482  //
483 
484  // Mark edges on intersection
485  const labelList& int1Edges = cutSurf1.intersectionEdges();
486 
487  boolList isIntersectionEdge1(cutSurf1.nEdges(), false);
488  forAll(int1Edges, intEdgeI)
489  {
490  label edgeI = int1Edges[intEdgeI];
491  isIntersectionEdge1[edgeI] = true;
492  }
493 
494  labelList faceZone1;
495  cutSurf1.markZones(isIntersectionEdge1, faceZone1);
496 
497 
498  // Check whether at least one of sides of intersection has been marked.
499  checkIncluded(cutSurf1, faceZone1, cutSurf1FaceI);
500 
501  // Subset zone which includes cutSurf2FaceI
502  boolList includedFaces1(cutSurf1.size(), false);
503 
504  forAll(faceZone1, faceI)
505  {
506  if (faceZone1[faceI] == faceZone1[cutSurf1FaceI])
507  {
508  includedFaces1[faceI] = true;
509  }
510  }
511 
512  // Subset to include only interesting part
513  labelList pointMap1;
514  labelList faceMap1;
515 
516  triSurface subSurf1
517  (
518  cutSurf1.subsetMesh
519  (
520  includedFaces1,
521  pointMap1,
522  faceMap1
523  )
524  );
525 
526 
527  //
528  // Mark faces of cutSurf2 that need to be kept by walking from includeFace2
529  // without crossing any edges of the intersection.
530  //
531 
532  // Mark edges and points on intersection
533  const labelList& int2Edges = cutSurf2.intersectionEdges();
534 
535  boolList isIntersectionEdge2(cutSurf2.nEdges(), false);
536  forAll(int2Edges, intEdgeI)
537  {
538  label edgeI = int2Edges[intEdgeI];
539  isIntersectionEdge2[edgeI] = true;
540  }
541 
542  labelList faceZone2;
543  cutSurf2.markZones(isIntersectionEdge2, faceZone2);
544 
545 
546  // Check whether at least one of sides of intersection has been marked.
547  checkIncluded(cutSurf2, faceZone2, cutSurf2FaceI);
548 
549  // Subset zone which includes cutSurf2FaceI
550  boolList includedFaces2(cutSurf2.size(), false);
551 
552  forAll(faceZone2, faceI)
553  {
554  if (faceZone2[faceI] == faceZone2[cutSurf2FaceI])
555  {
556  includedFaces2[faceI] = true;
557  }
558  }
559 
560  labelList pointMap2;
561  labelList faceMap2;
562 
563  triSurface subSurf2
564  (
565  cutSurf2.subsetMesh
566  (
567  includedFaces2,
568  pointMap2,
569  faceMap2
570  )
571  );
572 
573 
574  //
575  // Now match up the corresponding points on the intersection. The
576  // intersectedSurfaces will have the points resulting from the
577  // intersection last in their points and in the same
578  // order so we can use the pointMaps from the subsets to find them.
579  //
580  // We keep the vertices on the first surface and renumber those on the
581  // second one.
582 
583 
584  //
585  // points
586  //
587  pointField combinedPoints
588  (
589  subSurf1.nPoints()
590  + subSurf2.nPoints()
591  - (cutSurf2.nPoints() - cutSurf2.nSurfacePoints())
592  );
593 
594  // Copy points from subSurf1 and remember the labels of the ones in
595  // the intersection
596  labelList intersectionLabels
597  (
598  cutSurf1.nPoints() - cutSurf1.nSurfacePoints()
599  );
600 
601  label combinedPointI = 0;
602 
603  forAll(subSurf1.points(), pointI)
604  {
605  // Label in cutSurf
606  label cutSurfPointI = pointMap1[pointI];
607 
608  if (!cutSurf1.isSurfacePoint(cutSurfPointI))
609  {
610  // Label in original intersection is equal to the cutSurfPointI
611 
612  // Remember label in combinedPoints for intersection point.
613  intersectionLabels[cutSurfPointI] = combinedPointI;
614  }
615 
616  // Copy point
617  combinedPoints[combinedPointI++] = subSurf1.points()[pointI];
618  }
619 
620  // Append points from subSurf2 (if they are not intersection points)
621  // and construct mapping
622  labelList pointMap(subSurf2.nPoints());
623 
624  forAll(subSurf2.points(), pointI)
625  {
626  // Label in cutSurf
627  label cutSurfPointI = pointMap2[pointI];
628 
629  if (!cutSurf2.isSurfacePoint(cutSurfPointI))
630  {
631  // Lookup its label in combined point list.
632  pointMap[pointI] = intersectionLabels[cutSurfPointI];
633  }
634  else
635  {
636  pointMap[pointI] = combinedPointI;
637 
638  combinedPoints[combinedPointI++] = subSurf2.points()[pointI];
639  }
640  }
641 
642 
643  //
644  // patches
645  //
646 
647  labelList patchMap2;
648 
649  geometricSurfacePatchList combinedPatches
650  (
651  mergePatches
652  (
653  surf1,
654  surf2,
655  patchMap2
656  )
657  );
658 
659 
660  //
661  // faces
662  //
663 
664  List<labelledTri> combinedFaces(subSurf1.size() + subSurf2.size());
665 
666  faceMap_.setSize(combinedFaces.size());
667 
668  // Copy faces from subSurf1. No need for renumbering.
669  label combinedFaceI = 0;
670  forAll(subSurf1, faceI)
671  {
672  faceMap_[combinedFaceI] = faceMap1[faceI];
673  combinedFaces[combinedFaceI++] = subSurf1[faceI];
674  }
675 
676  // Copy and renumber faces from subSurf2.
677  forAll(subSurf2, faceI)
678  {
679  const labelledTri& f = subSurf2[faceI];
680 
681  faceMap_[combinedFaceI] = -faceMap2[faceI]-1;
682 
683  combinedFaces[combinedFaceI++] =
685  (
686  pointMap[f[0]],
687  pointMap[f[1]],
688  pointMap[f[2]],
689  patchMap2[f.region()]
690  );
691  }
692 
693  triSurface::operator=
694  (
695  triSurface
696  (
697  combinedFaces,
698  combinedPatches,
699  combinedPoints
700  )
701  );
702 }
703 
704 
705 // Construct from surfaces and boolean operation
707 (
708  const triSurface& surf1,
709  const triSurface& surf2,
710  const surfaceIntersection& inter,
711  const label booleanOp
712 )
713 :
714  triSurface(),
715  faceMap_()
716 {
717  if (debug)
718  {
719  Pout<< "booleanSurface : Testing surf1 and surf2" << endl;
720 
721  {
722  const labelListList& edgeFaces = surf1.edgeFaces();
723 
724  forAll(edgeFaces, edgeI)
725  {
726  const labelList& eFaces = edgeFaces[edgeI];
727 
728  if (eFaces.size() == 1)
729  {
731  << "surf1 is open surface at edge " << edgeI
732  << " verts:" << surf1.edges()[edgeI]
733  << " connected to faces " << eFaces << endl;
734  }
735  }
736  }
737  {
738  const labelListList& edgeFaces = surf2.edgeFaces();
739 
740  forAll(edgeFaces, edgeI)
741  {
742  const labelList& eFaces = edgeFaces[edgeI];
743 
744  if (eFaces.size() == 1)
745  {
747  << "surf2 is open surface at edge " << edgeI
748  << " verts:" << surf2.edges()[edgeI]
749  << " connected to faces " << eFaces << endl;
750  }
751  }
752  }
753  }
754 
755 
756  //
757  // Surface 1
758  //
759 
760  if (debug)
761  {
762  Pout<< "booleanSurface : Generating intersected surface for surf1"
763  << endl;
764  }
765 
766  // Add intersection to surface1 (retriangulates cut faces)
767  intersectedSurface cutSurf1(surf1, true, inter);
768 
769  if (debug)
770  {
771  Pout<< "booleanSurface : Generated cutSurf1: " << endl;
772  cutSurf1.writeStats(Pout);
773 
774  Pout<< "Writing to file cutSurf1.obj" << endl;
775  cutSurf1.write("cutSurf1.obj");
776  }
777 
778 
779  //
780  // Surface 2
781  //
782 
783  if (debug)
784  {
785  Pout<< "booleanSurface : Generating intersected surface for surf2"
786  << endl;
787  }
788 
789  // Add intersection to surface2
790  intersectedSurface cutSurf2(surf2, false, inter);
791 
792  if (debug)
793  {
794  Pout<< "booleanSurface : Generated cutSurf2: " << endl;
795  cutSurf2.writeStats(Pout);
796 
797  Pout<< "Writing to file cutSurf2.obj" << endl;
798  cutSurf2.write("cutSurf2.obj");
799  }
800 
801 
802  //
803  // patches
804  //
805 
806  labelList patchMap2;
807 
808  geometricSurfacePatchList combinedPatches
809  (
810  mergePatches
811  (
812  surf1,
813  surf2,
814  patchMap2
815  )
816  );
817 
818 
819  //
820  // Now match up the corresponding points on the intersection. The
821  // intersectedSurfaces will have the points resulting from the
822  // intersection first in their points and in the same
823  // order
824  //
825  // We keep the vertices on the first surface and renumber those on the
826  // second one.
827 
828  pointField combinedPoints(cutSurf1.nPoints() + cutSurf2.nSurfacePoints());
829 
830  // Copy all points from 1 and non-intersection ones from 2.
831 
832  label combinedPointI = 0;
833 
834  forAll(cutSurf1.points(), pointI)
835  {
836  combinedPoints[combinedPointI++] = cutSurf1.points()[pointI];
837  }
838 
839  for
840  (
841  label pointI = 0;
842  pointI < cutSurf2.nSurfacePoints();
843  pointI++
844  )
845  {
846  combinedPoints[combinedPointI++] = cutSurf2.points()[pointI];
847  }
848 
849  // Point order is now
850  // - 0.. cutSurf1.nSurfacePoints : original surf1 points
851  // - .. cutSurf1.nPoints : intersection points
852  // - .. cutSurf2.nSurfacePoints : original surf2 points
853 
854  if (debug)
855  {
856  Pout<< "booleanSurface : generated points:" << nl
857  << " 0 .. " << cutSurf1.nSurfacePoints()-1
858  << " : original surface1"
859  << nl
860  << " " << cutSurf1.nSurfacePoints()
861  << " .. " << cutSurf1.nPoints()-1
862  << " : intersection points"
863  << nl
864  << " " << cutSurf1.nPoints() << " .. "
865  << cutSurf2.nSurfacePoints()-1
866  << " : surface2 points"
867  << nl
868  << endl;
869  }
870 
871  // Copy faces. Faces from surface 1 keep vertex numbering and region info.
872  // Faces from 2 get vertices and region renumbered.
873  List<labelledTri> combinedFaces(cutSurf1.size() + cutSurf2.size());
874 
875  label combinedFaceI = 0;
876 
877  forAll(cutSurf1, faceI)
878  {
879  combinedFaces[combinedFaceI++] = cutSurf1[faceI];
880  }
881 
882  forAll(cutSurf2, faceI)
883  {
884  labelledTri& combinedTri = combinedFaces[combinedFaceI++];
885 
886  const labelledTri& tri = cutSurf2[faceI];
887 
888  forAll(tri, fp)
889  {
890  if (cutSurf2.isSurfacePoint(tri[fp]))
891  {
892  // Renumber. Surface2 points are after ones from surf 1.
893  combinedTri[fp] = tri[fp] + cutSurf1.nPoints();
894  }
895  else
896  {
897  // Is intersection.
898  combinedTri[fp] =
899  tri[fp]
900  - cutSurf2.nSurfacePoints()
901  + cutSurf1.nSurfacePoints();
902  }
903  }
904  combinedTri.region() = patchMap2[tri.region()];
905  }
906 
907 
908  // Now we have surface in combinedFaces and combinedPoints. Use
909  // booleanOp to determine which part of what to keep.
910 
911  // Construct addressing for whole part.
912  triSurface combinedSurf
913  (
914  combinedFaces,
915  combinedPatches,
916  combinedPoints
917  );
918 
919  if (debug)
920  {
921  Pout<< "booleanSurface : Generated combinedSurf: " << endl;
922  combinedSurf.writeStats(Pout);
923 
924  Pout<< "Writing to file combinedSurf.obj" << endl;
925  combinedSurf.write("combinedSurf.obj");
926  }
927 
928 
929  if (booleanOp == booleanSurface::ALL)
930  {
931  // Special case: leave surface multiply connected
932 
933  faceMap_.setSize(combinedSurf.size());
934 
935  label combinedFaceI = 0;
936 
937  forAll(cutSurf1, faceI)
938  {
939  faceMap_[combinedFaceI++] = cutSurf1.faceMap()[faceI];
940  }
941  forAll(cutSurf2, faceI)
942  {
943  faceMap_[combinedFaceI++] = -cutSurf2.faceMap()[faceI] - 1;
944  }
945 
946  triSurface::operator=(combinedSurf);
947 
948  return;
949  }
950 
951 
952  // Get outside point.
953  point outsidePoint = 2 * treeBoundBox(combinedSurf.localPoints()).span();
954 
955  //
956  // Linear search for nearest point on surface.
957  //
958 
959  const pointField& pts = combinedSurf.points();
960 
961  label minFaceI = -1;
962  pointHit minHit(false, vector::zero, GREAT, true);
963 
964  forAll(combinedSurf, faceI)
965  {
966  pointHit curHit = combinedSurf[faceI].nearestPoint(outsidePoint, pts);
967 
968  if (curHit.distance() < minHit.distance())
969  {
970  minHit = curHit;
971  minFaceI = faceI;
972  }
973  }
974 
975  if (debug)
976  {
977  Pout<< "booleanSurface : found for point:" << outsidePoint
978  << " nearest face:" << minFaceI
979  << " nearest point:" << minHit.rawPoint()
980  << endl;
981  }
982 
983  // Visibility/side of face:
984  // UNVISITED: unvisited
985  // OUTSIDE: visible from outside
986  // INSIDE: invisible from outside
987  labelList side(combinedSurf.size(), UNVISITED);
988 
989  // Walk face-edge-face and propagate inside/outside status.
990  propagateSide(combinedSurf, OUTSIDE, minFaceI, side);
991 
992 
993  // Depending on operation include certain faces.
994  // INTERSECTION: faces on inside of 1 and of 2
995  // UNION: faces on outside of 1 and of 2
996  // DIFFERENCE: faces on outside of 1 and inside of 2
997 
998  boolList include(combinedSurf.size(), false);
999 
1000  forAll(side, faceI)
1001  {
1002  if (side[faceI] == UNVISITED)
1003  {
1005  << "Face " << faceI << " has not been reached by walking from"
1006  << " nearest point " << minHit.rawPoint()
1007  << " nearest face " << minFaceI << exit(FatalError);
1008  }
1009  else if (side[faceI] == OUTSIDE)
1010  {
1011  if (booleanOp == booleanSurface::UNION)
1012  {
1013  include[faceI] = true;
1014  }
1015  else if (booleanOp == booleanSurface::INTERSECTION)
1016  {
1017  include[faceI] = false;
1018  }
1019  else // difference
1020  {
1021  include[faceI] = (faceI < cutSurf1.size()); // face from surf1
1022  }
1023  }
1024  else // inside
1025  {
1026  if (booleanOp == booleanSurface::UNION)
1027  {
1028  include[faceI] = false;
1029  }
1030  else if (booleanOp == booleanSurface::INTERSECTION)
1031  {
1032  include[faceI] = true;
1033  }
1034  else // difference
1035  {
1036  include[faceI] = (faceI >= cutSurf1.size()); // face from surf2
1037  }
1038  }
1039  }
1040 
1041  // Create subsetted surface
1042  labelList subToCombinedPoint;
1043  labelList subToCombinedFace;
1044  triSurface subSurf
1045  (
1046  combinedSurf.subsetMesh
1047  (
1048  include,
1049  subToCombinedPoint,
1050  subToCombinedFace
1051  )
1052  );
1053 
1054  // Create face map
1055  faceMap_.setSize(subSurf.size());
1056 
1057  forAll(subToCombinedFace, faceI)
1058  {
1059  // Get label in combinedSurf
1060  label combinedFaceI = subToCombinedFace[faceI];
1061 
1062  // First faces in combinedSurf come from cutSurf1.
1063 
1064  if (combinedFaceI < cutSurf1.size())
1065  {
1066  label cutSurf1Face = combinedFaceI;
1067 
1068  faceMap_[faceI] = cutSurf1.faceMap()[cutSurf1Face];
1069  }
1070  else
1071  {
1072  label cutSurf2Face = combinedFaceI - cutSurf1.size();
1073 
1074  faceMap_[faceI] = - cutSurf2.faceMap()[cutSurf2Face] - 1;
1075  }
1076  }
1077 
1078  // Orient outwards
1079  orientedSurface outSurf(subSurf);
1080 
1081  // Assign.
1082  triSurface::operator=(outSurf);
1083 }
1084 
1085 
1086 // ************************************************************************* //
meshTools.H
Foam::booleanSurface::propagateEdgeSide
static void propagateEdgeSide(const triSurface &surf, const label prevVert0, const label prevFaceI, const label prevState, const label edgeI, labelList &side)
On edgeI, coming from face prevFace, determines visibility/side of.
Definition: booleanSurface.C:197
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Return edge-face addressing.
Definition: PrimitivePatchTemplate.C:292
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
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatchTemplate.C:212
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::PointHit
This class describes the interaction of a face and a point. It carries the info of a successful hit a...
Definition: PointHit.H:51
Foam::orientedSurface
Given point flip all faces such that normals point in same direction.
Definition: orientedSurface.H:52
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
Foam::PointHit::rawPoint
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
Foam::PrimitivePatch::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatchTemplate.H:299
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::intersectedSurface::nSurfacePoints
label nSurfacePoints() const
Number of points from original surface.
Definition: intersectedSurface.H:280
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::booleanSurface::UNION
@ UNION
Definition: booleanSurface.H:165
Foam::PrimitivePatch::faceEdges
const labelListList & faceEdges() const
Return face-edge addressing.
Definition: PrimitivePatchTemplate.C:312
Foam::booleanSurface::ALL
@ ALL
Definition: booleanSurface.H:168
Foam::triSurface::patches
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:301
Foam::booleanSurface::mergePatches
static geometricSurfacePatchList mergePatches(const triSurface &surf1, const triSurface &surf2, labelList &patchMap2)
Generate combined patchList (returned). Sets patchMap to map from.
Definition: booleanSurface.C:141
Foam::booleanSurface::booleanOpType
booleanOpType
Enumeration listing the possible volume operator types.
Definition: booleanSurface.H:163
OFstream.H
Foam::boundBox::span
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:84
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::PointHit::distance
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
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::meshTools::findEdge
label findEdge(const edgeList &edges, const labelList &candidates, const label v0, const label v1)
Return edge among candidates that uses the two vertices.
Definition: meshTools.C:336
Foam::intersectedSurface
Given triSurface and intersection creates the intersected (properly triangulated) surface....
Definition: intersectedSurface.H:81
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::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
treeBoundBox.H
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::booleanSurface::index
static label index(const labelList &elems, const label elem)
Get label in elems of elem.
Definition: booleanSurface.C:100
Foam::triSurface::operator=
void operator=(const triSurface &)
Definition: triSurface.C:1122
Foam::intersectedSurface::faceMap
const labelList & faceMap() const
New to old.
Definition: intersectedSurface.H:274
Foam::PrimitivePatch::nPoints
label nPoints() const
Return number of points supporting patch faces.
Definition: PrimitivePatchTemplate.H:293
Foam::FatalError
error FatalError
intersectedSurface.H
Foam::intersectedSurface::intersectionEdges
const labelList & intersectionEdges() const
Labels of edges in *this which originate from 'cuts'.
Definition: intersectedSurface.H:268
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::triSurface::markZones
label markZones(const boolList &borderEdge, labelList &faceZone) const
(size and) fills faceZone with zone of face. Zone is area
Definition: triSurface.C:906
Foam::booleanSurface::propagateSide
static void propagateSide(const triSurface &surf, const label prevState, const label faceI, labelList &side)
Given in/outside status of face determines status for all.
Definition: booleanSurface.C:324
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
Foam::labelledTri::region
label region() const
Return region label.
Definition: labelledTriI.H:68
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
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::booleanSurface::checkIncluded
static void checkIncluded(const intersectedSurface &surf, const labelList &faceZone, const label includedFace)
Check whether subset of faces (from markZones) reaches up to.
Definition: booleanSurface.C:63
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::booleanSurface::booleanSurface
booleanSurface()
Construct null.
Definition: booleanSurface.C:388
Foam::labelledTri
Triangle with additional region number.
Definition: labelledTri.H:49
triSurfaceSearch.H
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::PrimitivePatch::localFaces
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatchTemplate.C:372
Foam::triSurface::sortedEdgeFaces
const labelListList & sortedEdgeFaces() const
Return edge-face addressing sorted (for edges with more than.
Definition: triSurface.C:766
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::intersectedSurface::isSurfacePoint
bool isSurfacePoint(const label pointI) const
Is point coming from original surface?
Definition: intersectedSurface.H:286
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
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
orientedSurface.H
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
Foam::booleanSurface::findEdge
static label findEdge(const edgeList &edges, const labelList &edgeLabels, const edge &e)
Find index of edge e in subset edgeLabels.
Definition: booleanSurface.C:117