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-2013 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 \*---------------------------------------------------------------------------*/
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 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 // Check whether at least one of faces connected to the intersection has been
45 // marked.
46 void Foam::booleanSurface::checkIncluded
47 (
48  const intersectedSurface& surf,
49  const labelList& faceZone,
50  const label includedFace
51 )
52 {
53  forAll(surf.intersectionEdges(), intEdgeI)
54  {
55  label edgeI = surf.intersectionEdges()[intEdgeI];
56 
57  const labelList& myFaces = surf.edgeFaces()[edgeI];
58 
59  bool usesIncluded = false;
60 
61  forAll(myFaces, myFaceI)
62  {
63  if (faceZone[myFaces[myFaceI]] == faceZone[includedFace])
64  {
65  usesIncluded = true;
66 
67  break;
68  }
69  }
70 
71  if (!usesIncluded)
72  {
74  (
75  "booleanSurface::checkIncluded(const intersectedSurface&"
76  ", const labelList&, const label)"
77  ) << "None of the faces reachable from face " << includedFace
78  << " connects to the intersection."
79  << exit(FatalError);
80  }
81  }
82 }
83 
84 
85 // Linear lookup
86 Foam::label Foam::booleanSurface::index
87 (
88  const labelList& elems,
89  const label elem
90 )
91 {
92  forAll(elems, elemI)
93  {
94  if (elems[elemI] == elem)
95  {
96  return elemI;
97  }
98  }
99  return -1;
100 }
101 
102 
103 Foam::label Foam::booleanSurface::findEdge
104 (
105  const edgeList& edges,
106  const labelList& edgeLabels,
107  const edge& e
108 )
109 {
110  forAll(edgeLabels, edgeLabelI)
111  {
112  if (edges[edgeLabels[edgeLabelI]] == e)
113  {
114  return edgeLabels[edgeLabelI];
115  }
116  }
118  (
119  "booleanSurface::findEdge(const edgeList&, const labelList&"
120  ", const edge&)"
121  ) << "Cannot find edge " << e << " in edges " << edgeLabels
122  << abort(FatalError);
123 
124  return -1;
125 }
126 
127 
128 // Generate combined patchList (returned). Sets patchMap to map from surf2
129 // region numbers into combined patchList
130 Foam::geometricSurfacePatchList Foam::booleanSurface::mergePatches
131 (
132  const triSurface& surf1,
133  const triSurface& surf2,
134  labelList& patchMap2
135 )
136 {
137  // Size too big.
138  geometricSurfacePatchList combinedPatches
139  (
140  surf1.patches().size()
141  + surf2.patches().size()
142  );
143 
144  // Copy all patches of surf1
145  label combinedPatchI = 0;
146  forAll(surf1.patches(), patchI)
147  {
148  combinedPatches[combinedPatchI++] = surf1.patches()[patchI];
149  }
150 
151  // (inefficiently) add unique patches from surf2
152  patchMap2.setSize(surf2.patches().size());
153 
154  forAll(surf2.patches(), patch2I)
155  {
156  label index = -1;
157 
158  forAll(surf1.patches(), patch1I)
159  {
160  if (surf1.patches()[patch1I] == surf2.patches()[patch2I])
161  {
162  index = patch1I;
163 
164  break;
165  }
166  }
167 
168  if (index == -1)
169  {
170  combinedPatches[combinedPatchI] = surf2.patches()[patch2I];
171  patchMap2[patch2I] = combinedPatchI;
172  combinedPatchI++;
173  }
174  else
175  {
176  patchMap2[patch2I] = index;
177  }
178  }
179 
180  combinedPatches.setSize(combinedPatchI);
181 
182  return combinedPatches;
183 }
184 
185 
186 void Foam::booleanSurface::propagateEdgeSide
187 (
188  const triSurface& surf,
189  const label prevVert0,
190  const label prevFaceI,
191  const label prevState,
192  const label edgeI,
193  labelList& side
194 )
195 {
196  const labelList& eFaces = surf.sortedEdgeFaces()[edgeI];
197 
198  // Simple case. Propagate side.
199  if (eFaces.size() == 2)
200  {
201  forAll(eFaces, eFaceI)
202  {
203  propagateSide
204  (
205  surf,
206  prevState,
207  eFaces[eFaceI],
208  side
209  );
210  }
211  }
212 
213 
214  if (((eFaces.size() % 2) == 1) && (eFaces.size() != 1))
215  {
217  (
218  "booleanSurface::propagateEdgeSide(const triSurface&,"
219  "const label, const label, const label, const label,"
220  " labelList&)"
221  ) << "Don't know how to handle edges with odd number of faces"
222  << endl
223  << "edge:" << edgeI << " vertices:" << surf.edges()[edgeI]
224  << " coming from face:" << prevFaceI
225  << " edgeFaces:" << eFaces << abort(FatalError);
226  }
227 
228 
229  // Get position of face in edgeFaces
230  label ind = index(eFaces, prevFaceI);
231 
232  // Determine orientation of faces around edge prevVert0
233  // (might be opposite of edge)
234  const edge& e = surf.edges()[edgeI];
235 
236  // Get next face to include
237  label nextInd;
238  label prevInd;
239 
240  if (e.start() == prevVert0)
241  {
242  // Edge (and hence eFaces) in same order as prevVert0.
243  // Take next face from sorted list
244  nextInd = eFaces.fcIndex(ind);
245  prevInd = eFaces.rcIndex(ind);
246  }
247  else
248  {
249  // Take previous face from sorted neighbours
250  nextInd = eFaces.rcIndex(ind);
251  prevInd = eFaces.fcIndex(ind);
252  }
253 
254 
255  if (prevState == OUTSIDE)
256  {
257  // Coming from outside. nextInd is outside, rest is inside.
258 
259  forAll(eFaces, eFaceI)
260  {
261  if (eFaceI != ind)
262  {
263  label nextState;
264 
265  if (eFaceI == nextInd)
266  {
267  nextState = OUTSIDE;
268  }
269  else
270  {
271  nextState = INSIDE;
272  }
273 
274  propagateSide
275  (
276  surf,
277  nextState,
278  eFaces[eFaceI],
279  side
280  );
281  }
282  }
283  }
284  else
285  {
286  // Coming from inside. prevInd is inside as well, rest is outside.
287 
288  forAll(eFaces, eFaceI)
289  {
290  if (eFaceI != ind)
291  {
292  label nextState;
293 
294  if (eFaceI == prevInd)
295  {
296  nextState = INSIDE;
297  }
298  else
299  {
300  nextState = OUTSIDE;
301  }
302 
303  propagateSide
304  (
305  surf,
306  nextState,
307  eFaces[eFaceI],
308  side
309  );
310  }
311  }
312  }
313 }
314 
315 
316 // Face-edge walk. Determines inside/outside for all faces connected to an edge.
317 void Foam::booleanSurface::propagateSide
318 (
319  const triSurface& surf,
320  const label prevState,
321  const label faceI,
322  labelList& side
323 )
324 {
325  if (side[faceI] == UNVISITED)
326  {
327  side[faceI] = prevState;
328 
329  const labelledTri& tri = surf.localFaces()[faceI];
330 
331  // Get copy of face labels
332  label a = tri[0];
333  label b = tri[1];
334  label c = tri[2];
335 
336  // Go and visit my edges' face-neighbours.
337 
338  const labelList& myEdges = surf.faceEdges()[faceI];
339 
340  label edgeAB = findEdge(surf.edges(), myEdges, edge(a, b));
341 
342  propagateEdgeSide
343  (
344  surf,
345  a,
346  faceI,
347  prevState,
348  edgeAB,
349  side
350  );
351 
352  label edgeBC = findEdge(surf.edges(), myEdges, edge(b, c));
353 
354  propagateEdgeSide
355  (
356  surf,
357  b,
358  faceI,
359  prevState,
360  edgeBC,
361  side
362  );
363 
364  label edgeCA = findEdge(surf.edges(), myEdges, edge(c, a));
365 
366  propagateEdgeSide
367  (
368  surf,
369  c,
370  faceI,
371  prevState,
372  edgeCA,
373  side
374  );
375  }
376 }
377 
378 
379 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
380 
381 // Null constructor
383 :
384  triSurface()
385 {}
386 
387 
388 // Construct from surfaces and face to include for every surface
390 (
391  const triSurface& surf1,
392  const triSurface& surf2,
393  const surfaceIntersection& inter,
394  const label includeFace1,
395  const label includeFace2
396 )
397 :
398  triSurface(),
399  faceMap_()
400 {
401  if (debug)
402  {
403  Pout<< "booleanSurface : Generating intersected surface for surf1"
404  << endl;
405  }
406 
407  // Add intersection to surface1 (retriangulates cut faces)
408  intersectedSurface cutSurf1(surf1, true, inter);
409 
410 
411  if (debug)
412  {
413  Pout<< "booleanSurface : Generated cutSurf1: " << endl;
414  cutSurf1.writeStats(Pout);
415 
416  Pout<< "Writing to file cutSurf1.obj" << endl;
417  cutSurf1.write("cutSurf1.obj");
418  }
419 
420  if (debug)
421  {
422  Pout<< "booleanSurface : Generating intersected surface for surf2"
423  << endl;
424  }
425 
426  // Add intersection to surface2
427  intersectedSurface cutSurf2(surf2, false, inter);
428 
429  if (debug)
430  {
431  Pout<< "booleanSurface : Generated cutSurf2: " << endl;
432  cutSurf2.writeStats(Pout);
433 
434  Pout<< "Writing to file cutSurf2.obj" << endl;
435  cutSurf2.write("cutSurf2.obj");
436  }
437 
438 
439  // Find (first) face of cutSurf1 that originates from includeFace1
440  label cutSurf1FaceI = index(cutSurf1.faceMap(), includeFace1);
441 
442  if (debug)
443  {
444  Pout<< "cutSurf1 : starting to fill from face:" << cutSurf1FaceI
445  << endl;
446  }
447 
448  if (cutSurf1FaceI == -1)
449  {
451  (
452  "booleanSurface(const triSurfaceSearch&"
453  ", const label, const triSurfaceSearch&, const label)"
454  ) << "Did not find face with label " << includeFace1
455  << " in intersectedSurface."
456  << exit(FatalError);
457  }
458 
459  // Find (first) face of cutSurf2 that originates from includeFace1
460  label cutSurf2FaceI = index(cutSurf2.faceMap(), includeFace2);
461 
462  if (debug)
463  {
464  Pout<< "cutSurf2 : starting to fill from face:" << cutSurf2FaceI
465  << endl;
466  }
467  if (cutSurf2FaceI == -1)
468  {
470  (
471  "booleanSurface(const triSurfaceSearch&"
472  ", const label, const triSurfaceSearch&, const label)"
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 
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  {
730  WarningIn("booleanSurface::booleanSurface")
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  {
746  WarningIn("booleanSurface::booleanSurface")
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  {
1004  FatalErrorIn
1005  (
1006  "booleanSurface::booleanSurface"
1007  "(const triSurfaceSearch&, const triSurfaceSearch&"
1008  ", const label booleanOp)"
1009  ) << "Face " << faceI << " has not been reached by walking from"
1010  << " nearest point " << minHit.rawPoint()
1011  << " nearest face " << minFaceI << exit(FatalError);
1012  }
1013  else if (side[faceI] == OUTSIDE)
1014  {
1015  if (booleanOp == booleanSurface::UNION)
1016  {
1017  include[faceI] = true;
1018  }
1019  else if (booleanOp == booleanSurface::INTERSECTION)
1020  {
1021  include[faceI] = false;
1022  }
1023  else // difference
1024  {
1025  include[faceI] = (faceI < cutSurf1.size()); // face from surf1
1026  }
1027  }
1028  else // inside
1029  {
1030  if (booleanOp == booleanSurface::UNION)
1031  {
1032  include[faceI] = false;
1033  }
1034  else if (booleanOp == booleanSurface::INTERSECTION)
1035  {
1036  include[faceI] = true;
1037  }
1038  else // difference
1039  {
1040  include[faceI] = (faceI >= cutSurf1.size()); // face from surf2
1041  }
1042  }
1043  }
1044 
1045  // Create subsetted surface
1046  labelList subToCombinedPoint;
1047  labelList subToCombinedFace;
1048  triSurface subSurf
1049  (
1050  combinedSurf.subsetMesh
1051  (
1052  include,
1053  subToCombinedPoint,
1054  subToCombinedFace
1055  )
1056  );
1057 
1058  // Create face map
1059  faceMap_.setSize(subSurf.size());
1060 
1061  forAll(subToCombinedFace, faceI)
1062  {
1063  // Get label in combinedSurf
1064  label combinedFaceI = subToCombinedFace[faceI];
1065 
1066  // First faces in combinedSurf come from cutSurf1.
1067 
1068  if (combinedFaceI < cutSurf1.size())
1069  {
1070  label cutSurf1Face = combinedFaceI;
1071 
1072  faceMap_[faceI] = cutSurf1.faceMap()[cutSurf1Face];
1073  }
1074  else
1075  {
1076  label cutSurf2Face = combinedFaceI - cutSurf1.size();
1077 
1078  faceMap_[faceI] = - cutSurf2.faceMap()[cutSurf2Face] - 1;
1079  }
1080  }
1081 
1082  // Orient outwards
1083  orientedSurface outSurf(subSurf);
1084 
1085  // Assign.
1086  triSurface::operator=(outSurf);
1087 }
1088 
1089 
1090 // ************************************************************************* //
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
bool isSurfacePoint(const label pointI) const
Is point coming from original surface?
labelList f(nPoints)
Basic surface-surface intersection description. Constructed from two surfaces it creates a descriptio...
const labelListList & edgeFaces() const
Return edge-face addressing.
Given point flip all faces such that normals point in same direction.
const labelList & faceMap() const
New to old.
Given triSurface and intersection creates the intersected (properly triangulated) surface...
Triangulated surface description with patch information.
Definition: triSurface.H:57
label markZones(const boolList &borderEdge, labelList &faceZone) const
(size and) fills faceZone with zone of face. Zone is area
Definition: triSurface.C:921
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
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
triSurface()
Construct null.
Definition: triSurface.C:611
friend Ostream & operator(Ostream &, const UList< T > &)
List< geometricSurfacePatch > geometricSurfacePatchList
label nPoints() const
Return number of points supporting patch faces.
Namespace for OpenFOAM.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Triangle with additional region number.
Definition: labelledTri.H:49
static const char nl
Definition: Ostream.H:260
label region() const
Return region label.
Definition: labelledTriI.H:68
void setSize(const label)
Reset size of List.
Definition: List.C:318
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
#define WarningIn(functionName)
Report a warning using Foam::Warning.
label nSurfacePoints() const
Number of points from original surface.
#define forAll(list, i)
Definition: UList.H:421
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
triSurface subsetMesh(const boolList &include, labelList &pointMap, labelList &faceMap) const
Return new surface. Returns pointMap, faceMap from.
Definition: triSurface.C:1016
const labelList & intersectionEdges() const
Labels of edges in *this which originate from &#39;cuts&#39;.
booleanSurface()
Construct null.
void operator=(const triSurface &)
Definition: triSurface.C:1125
errorManip< error > abort(error &err)
Definition: errorManip.H:131
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
List< edge > edgeList
Definition: edgeList.H:38
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
void write(Ostream &) const
Write to Ostream in simple FOAM format.
Definition: triSurface.C:1067
error FatalError
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const Vector zero
Definition: Vector.H:80
label nEdges() const
Return number of edges in patch.
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:343
void writeStats(Ostream &) const
Write some statistics.
Definition: triSurface.C:1092
const Field< PointType > & points() const
Return reference to global points.
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53