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-2016 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 
29  Generates the extendedFeatureEdgeMesh for the interface between a boolean
30  operation on two surfaces. Assumes that the orientation of the surfaces is
31  correct:
32 
33  - if the operation is union or intersection, that both surface's normals
34  (n) have the same orientation with respect to a point, i.e. surfaces and b
35  are orientated the same with respect to point x:
36 
37  @verbatim
38  _______
39  | |--> n
40  | ___|___ x
41  |a | | |--> n
42  |___|___| b|
43  | |
44  |_______|
45 
46  @endverbatim
47 
48  - if the operation is a subtraction, the surfaces should be oppositely
49  oriented with respect to a point, i.e. for (a - b), then b's orientation
50  should be such that x is "inside", and a's orientation such that x is
51  "outside"
52 
53  @verbatim
54  _______
55  | |--> n
56  | ___|___ x
57  |a | | |
58  |___|___| b|
59  | n <--|
60  |_______|
61 
62  @endverbatim
63 
64  When the operation is peformed - for union, all of the edges generates where
65  one surfaces cuts another are all "internal" for union, and "external" for
66  intersection, b - a and a - b. This has been assumed, formal (dis)proof is
67  invited.
68 
69 \*---------------------------------------------------------------------------*/
70 
71 #include "triSurface.H"
72 #include "argList.H"
73 #include "Time.H"
74 #include "featureEdgeMesh.H"
76 #include "triSurfaceSearch.H"
77 #include "OFstream.H"
78 #include "booleanSurface.H"
79 #include "edgeIntersections.H"
80 #include "meshTools.H"
81 #include "labelPair.H"
82 
83 using namespace Foam;
84 
85 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87 bool intersectSurfaces
88 (
89  triSurface& surf,
90  edgeIntersections& edgeCuts
91 )
92 {
93  bool hasMoved = false;
94 
95  for (label iter = 0; iter < 10; iter++)
96  {
97  Info<< "Determining intersections of surface edges with itself" << endl;
98 
99  // Determine surface edge intersections. Allow surface to be moved.
100 
101  // Number of iterations needed to resolve degenerates
102  label nIters = 0;
103  {
104  triSurfaceSearch querySurf(surf);
105 
106  scalarField surfPointTol
107  (
109  );
110 
111  // Determine raw intersections
112  edgeCuts = edgeIntersections
113  (
114  surf,
115  querySurf,
116  surfPointTol
117  );
118 
119  // Shuffle a bit to resolve degenerate edge-face hits
120  {
121  pointField points(surf.points());
122 
123  nIters =
124  edgeCuts.removeDegenerates
125  (
126  5, // max iterations
127  surf,
128  querySurf,
129  surfPointTol,
130  points // work array
131  );
132 
133  if (nIters != 0)
134  {
135  // Update geometric quantities
136  surf.movePoints(points);
137  hasMoved = true;
138  }
139  }
140  }
141  }
142 
143  if (hasMoved)
144  {
145  fileName newFile("surf.obj");
146  Info<< "Surface has been moved. Writing to " << newFile << endl;
147  surf.write(newFile);
148  }
149 
150  return hasMoved;
151 }
152 
153 
154 // Keep on shuffling surface points until no more degenerate intersections.
155 // Moves both surfaces and updates set of edge cuts.
156 bool intersectSurfaces
157 (
158  triSurface& surf1,
159  edgeIntersections& edgeCuts1,
160  triSurface& surf2,
161  edgeIntersections& edgeCuts2
162 )
163 {
164  bool hasMoved1 = false;
165  bool hasMoved2 = false;
166 
167  for (label iter = 0; iter < 10; iter++)
168  {
169  Info<< "Determining intersections of surf1 edges with surf2"
170  << " faces" << endl;
171 
172  // Determine surface1 edge intersections. Allow surface to be moved.
173 
174  // Number of iterations needed to resolve degenerates
175  label nIters1 = 0;
176  {
177  triSurfaceSearch querySurf2(surf2);
178 
179  scalarField surf1PointTol
180  (
182  );
183 
184  // Determine raw intersections
185  edgeCuts1 = edgeIntersections
186  (
187  surf1,
188  querySurf2,
189  surf1PointTol
190  );
191 
192  // Shuffle a bit to resolve degenerate edge-face hits
193  {
194  pointField points1(surf1.points());
195 
196  nIters1 =
197  edgeCuts1.removeDegenerates
198  (
199  5, // max iterations
200  surf1,
201  querySurf2,
202  surf1PointTol,
203  points1 // work array
204  );
205 
206  if (nIters1 != 0)
207  {
208  // Update geometric quantities
209  surf1.movePoints(points1);
210  hasMoved1 = true;
211  }
212  }
213  }
214 
215  Info<< "Determining intersections of surf2 edges with surf1"
216  << " faces" << endl;
217 
218  label nIters2 = 0;
219  {
220  triSurfaceSearch querySurf1(surf1);
221 
222  scalarField surf2PointTol
223  (
225  );
226 
227  // Determine raw intersections
228  edgeCuts2 = edgeIntersections
229  (
230  surf2,
231  querySurf1,
232  surf2PointTol
233  );
234 
235  // Shuffle a bit to resolve degenerate edge-face hits
236  {
237  pointField points2(surf2.points());
238 
239  nIters2 =
240  edgeCuts2.removeDegenerates
241  (
242  5, // max iterations
243  surf2,
244  querySurf1,
245  surf2PointTol,
246  points2 // work array
247  );
248 
249  if (nIters2 != 0)
250  {
251  // Update geometric quantities
252  surf2.movePoints(points2);
253  hasMoved2 = true;
254  }
255  }
256  }
257 
258  if (nIters1 == 0 && nIters2 == 0)
259  {
260  Info<< "** Resolved all intersections to be proper edge-face pierce"
261  << endl;
262  break;
263  }
264  }
265 
266  if (hasMoved1)
267  {
268  fileName newFile("surf1.obj");
269  Info<< "Surface 1 has been moved. Writing to " << newFile
270  << endl;
271  surf1.write(newFile);
272  }
273 
274  if (hasMoved2)
275  {
276  fileName newFile("surf2.obj");
277  Info<< "Surface 2 has been moved. Writing to " << newFile
278  << endl;
279  surf2.write(newFile);
280  }
281 
282  return hasMoved1 || hasMoved2;
283 }
284 
285 
286 label calcNormalDirection
287 (
288  const vector& normal,
289  const vector& otherNormal,
290  const vector& edgeDir,
291  const vector& faceCentre,
292  const vector& pointOnEdge
293 )
294 {
295  vector cross = (normal ^ edgeDir);
296  cross /= mag(cross);
297 
298  vector fC0tofE0 = faceCentre - pointOnEdge;
299  fC0tofE0 /= mag(fC0tofE0);
300 
301  label nDir = ((cross & fC0tofE0) > 0.0 ? 1 : -1);
302 
303  nDir *= ((otherNormal & fC0tofE0) > 0.0 ? -1 : 1);
304 
305  return nDir;
306 }
307 
308 
309 void calcEdgeCuts
310 (
311  triSurface& surf1,
312  triSurface& surf2,
313  const bool perturb,
314  edgeIntersections& edge1Cuts,
315  edgeIntersections& edge2Cuts
316 )
317 {
318  if (perturb)
319  {
320  intersectSurfaces
321  (
322  surf1,
323  edge1Cuts,
324  surf2,
325  edge2Cuts
326  );
327  }
328  else
329  {
330  triSurfaceSearch querySurf2(surf2);
331 
332  Info<< "Determining intersections of surf1 edges with surf2 faces"
333  << endl;
334 
335  edge1Cuts = edgeIntersections
336  (
337  surf1,
338  querySurf2,
340  );
341 
342  triSurfaceSearch querySurf1(surf1);
343 
344  Info<< "Determining intersections of surf2 edges with surf1 faces"
345  << endl;
346 
347  edge2Cuts = edgeIntersections
348  (
349  surf2,
350  querySurf1,
352  );
353  }
354 }
355 
356 
357 void calcFeaturePoints(const pointField& points, const edgeList& edges)
358 {
359  edgeMesh eMesh(points, edges);
360 
361  const labelListList& pointEdges = eMesh.pointEdges();
362 
363  // Get total number of feature points
364  label nFeaturePoints = 0;
365  forAll(pointEdges, pI)
366  {
367  const labelList& pEdges = pointEdges[pI];
368 
369  if (pEdges.size() == 1)
370  {
371  nFeaturePoints++;
372  }
373  }
374 
375 
376  // Calculate addressing from feature point to cut point and cut edge
377  labelList featurePointToCutPoint(nFeaturePoints);
378  labelList featurePointToCutEdge(nFeaturePoints);
379 
380  label nFeatPts = 0;
381  forAll(pointEdges, pI)
382  {
383  const labelList& pEdges = pointEdges[pI];
384 
385  if (pEdges.size() == 1)
386  {
387  featurePointToCutPoint[nFeatPts] = pI;
388  featurePointToCutEdge[nFeatPts] = pEdges[0];
389  nFeatPts++;
390  }
391  }
392 }
393 
394 
395 int main(int argc, char *argv[])
396 {
398  argList::validArgs.append("action");
399  argList::validArgs.append("surface file");
400  argList::validArgs.append("surface file");
401 
403  (
404  "surf1Baffle",
405  "Mark surface 1 as a baffle"
406  );
407 
409  (
410  "surf2Baffle",
411  "Mark surface 2 as a baffle"
412  );
413 
415  (
416  "perturb",
417  "Perturb surface points to escape degenerate intersections"
418  );
419 
421  (
422  "invertedSpace",
423  "do the surfaces have inverted space orientation, "
424  "i.e. a point at infinity is considered inside. "
425  "This is only sensible for union and intersection."
426  );
427 
428  #include "setRootCase.H"
429  #include "createTime.H"
430 
431  word action(args.args()[1]);
432 
434  validActions.insert("intersection", booleanSurface::INTERSECTION);
435  validActions.insert("union", booleanSurface::UNION);
436  validActions.insert("difference", booleanSurface::DIFFERENCE);
437 
438  if (!validActions.found(action))
439  {
441  << "Unsupported action " << action << endl
442  << "Supported actions:" << validActions.toc() << abort(FatalError);
443  }
444 
445  fileName surf1Name(args.args()[2]);
446  Info<< "Reading surface " << surf1Name << endl;
447  triSurface surf1(surf1Name);
448 
449  Info<< surf1Name << " statistics:" << endl;
450  surf1.writeStats(Info);
451  Info<< endl;
452 
453  fileName surf2Name(args[3]);
454  Info<< "Reading surface " << surf2Name << endl;
455  triSurface surf2(surf2Name);
456 
457  Info<< surf2Name << " statistics:" << endl;
458  surf2.writeStats(Info);
459  Info<< endl;
460 
461  const bool surf1Baffle = args.optionFound("surf1Baffle");
462  const bool surf2Baffle = args.optionFound("surf2Baffle");
463 
464  edgeIntersections edge1Cuts;
465  edgeIntersections edge2Cuts;
466 
467  bool invertedSpace = args.optionFound("invertedSpace");
468 
469  if (invertedSpace && validActions[action] == booleanSurface::DIFFERENCE)
470  {
472  << "Inverted space only makes sense for union or intersection."
473  << exit(FatalError);
474  }
475 
476  // Calculate the points where the edges are cut by the other surface
477  calcEdgeCuts
478  (
479  surf1,
480  surf2,
481  args.optionFound("perturb"),
482  edge1Cuts,
483  edge2Cuts
484  );
485 
486  // Determine intersection edges from the edge cuts
487  surfaceIntersection inter
488  (
489  surf1,
490  edge1Cuts,
491  surf2,
492  edge2Cuts
493  );
494 
495  fileName sFeatFileName
496  (
497  surf1Name.lessExt().name()
498  + "_"
499  + surf2Name.lessExt().name()
500  + "_"
501  + action
502  );
503 
504  label nFeatEds = inter.cutEdges().size();
505 
506  DynamicList<vector> normals(2*nFeatEds);
507  vectorField edgeDirections(nFeatEds, Zero);
509  (
510  2*nFeatEds
511  );
512  List<DynamicList<label>> edgeNormals(nFeatEds);
513  List<DynamicList<label>> normalDirections(nFeatEds);
514 
515  forAllConstIter(labelPairLookup, inter.facePairToEdge(), iter)
516  {
517  const label& cutEdgeI = iter();
518  const labelPair& facePair = iter.key();
519 
520  const edge& fE = inter.cutEdges()[cutEdgeI];
521 
522  const vector& norm1 = surf1.faceNormals()[facePair.first()];
523  const vector& norm2 = surf2.faceNormals()[facePair.second()];
524 
525  DynamicList<label>& eNormals = edgeNormals[cutEdgeI];
526  DynamicList<label>& nDirections = normalDirections[cutEdgeI];
527 
528  edgeDirections[cutEdgeI] = fE.vec(inter.cutPoints());
529 
530  normals.append(norm1);
531  eNormals.append(normals.size() - 1);
532 
533  if (surf1Baffle)
534  {
535  normalVolumeTypes.append(extendedFeatureEdgeMesh::BOTH);
536 
537  nDirections.append(1);
538  }
539  else
540  {
541  normalVolumeTypes.append(extendedFeatureEdgeMesh::INSIDE);
542  nDirections.append
543  (
544  calcNormalDirection
545  (
546  norm1,
547  norm2,
548  edgeDirections[cutEdgeI],
549  surf1[facePair.first()].centre(surf1.points()),
550  inter.cutPoints()[fE.start()]
551  )
552  );
553  }
554 
555  normals.append(norm2);
556  eNormals.append(normals.size() - 1);
557 
558  if (surf2Baffle)
559  {
560  normalVolumeTypes.append(extendedFeatureEdgeMesh::BOTH);
561 
562  nDirections.append(1);
563  }
564  else
565  {
566  normalVolumeTypes.append(extendedFeatureEdgeMesh::INSIDE);
567 
568  nDirections.append
569  (
570  calcNormalDirection
571  (
572  norm2,
573  norm1,
574  edgeDirections[cutEdgeI],
575  surf2[facePair.second()].centre(surf2.points()),
576  inter.cutPoints()[fE.start()]
577  )
578  );
579  }
580 
581 
582  if (surf1Baffle)
583  {
584  normals.append(norm2);
585 
586  if (surf2Baffle)
587  {
588  normalVolumeTypes.append(extendedFeatureEdgeMesh::BOTH);
589 
590  nDirections.append(1);
591  }
592  else
593  {
594  normalVolumeTypes.append(extendedFeatureEdgeMesh::INSIDE);
595 
596  nDirections.append
597  (
598  calcNormalDirection
599  (
600  norm2,
601  norm1,
602  edgeDirections[cutEdgeI],
603  surf2[facePair.second()].centre(surf2.points()),
604  inter.cutPoints()[fE.start()]
605  )
606  );
607  }
608 
609  eNormals.append(normals.size() - 1);
610  }
611 
612  if (surf2Baffle)
613  {
614  normals.append(norm1);
615 
616  if (surf1Baffle)
617  {
618  normalVolumeTypes.append(extendedFeatureEdgeMesh::BOTH);
619 
620  nDirections.append(1);
621  }
622  else
623  {
624  normalVolumeTypes.append(extendedFeatureEdgeMesh::INSIDE);
625 
626  nDirections.append
627  (
628  calcNormalDirection
629  (
630  norm1,
631  norm2,
632  edgeDirections[cutEdgeI],
633  surf1[facePair.first()].centre(surf1.points()),
634  inter.cutPoints()[fE.start()]
635  )
636  );
637  }
638 
639  eNormals.append(normals.size() - 1);
640  }
641  }
642 
643 
644  label internalStart = -1;
645  label nIntOrExt = 0;
646  label nFlat = 0;
647  label nOpen = 0;
648  label nMultiple = 0;
649 
650  forAll(edgeNormals, eI)
651  {
652  label nEdNorms = edgeNormals[eI].size();
653 
654  if (nEdNorms == 1)
655  {
656  nOpen++;
657  }
658  else if (nEdNorms == 2)
659  {
660  const vector& n0(normals[edgeNormals[eI][0]]);
661  const vector& n1(normals[edgeNormals[eI][1]]);
662 
664  {
665  nFlat++;
666  }
667  else
668  {
669  nIntOrExt++;
670  }
671  }
672  else if (nEdNorms > 2)
673  {
674  nMultiple++;
675  }
676  }
677 
678  if (validActions[action] == booleanSurface::UNION)
679  {
680  if (!invertedSpace)
681  {
682  // All edges are internal
683  internalStart = 0;
684  }
685  else
686  {
687  // All edges are external
688  internalStart = nIntOrExt;
689  }
690  }
691  else if (validActions[action] == booleanSurface::INTERSECTION)
692  {
693  if (!invertedSpace)
694  {
695  // All edges are external
696  internalStart = nIntOrExt;
697  }
698  else
699  {
700  // All edges are internal
701  internalStart = 0;
702  }
703  }
704  else if (validActions[action] == booleanSurface::DIFFERENCE)
705  {
706  // All edges are external
707  internalStart = nIntOrExt;
708  }
709  else
710  {
712  << "Unsupported booleanSurface:booleanOpType and space "
713  << action << " " << invertedSpace
714  << abort(FatalError);
715  }
716 
717  // There are no feature points supported by surfaceIntersection
718  // Flat, open or multiple edges are assumed to be impossible
719  // Region edges are not explicitly supported by surfaceIntersection
720 
721  vectorField normalsTmp(normals);
723  (
724  normalVolumeTypes
725  );
726  labelListList edgeNormalsTmp(edgeNormals.size());
727  forAll(edgeNormalsTmp, i)
728  {
729  edgeNormalsTmp[i] = edgeNormals[i];
730  }
731  labelListList normalDirectionsTmp(normalDirections.size());
732  forAll(normalDirectionsTmp, i)
733  {
734  normalDirectionsTmp[i] = normalDirections[i];
735  }
736 
737  calcFeaturePoints(inter.cutPoints(), inter.cutEdges());
738 
740  (
741  IOobject
742  (
743  sFeatFileName + ".extendedFeatureEdgeMesh",
744  runTime.constant(),
745  "extendedFeatureEdgeMesh",
746  runTime,
749  ),
750  inter.cutPoints(),
751  inter.cutEdges(),
752 
753  0, // concaveStart,
754  0, // mixedStart,
755  0, // nonFeatureStart,
756 
757  internalStart, // internalStart,
758  nIntOrExt, // flatStart,
759  nIntOrExt + nFlat, // openStart,
760  nIntOrExt + nFlat + nOpen, // multipleStart,
761 
762  normalsTmp,
763  normalVolumeTypesTmp,
764  edgeDirections,
765  normalDirectionsTmp,
766  edgeNormalsTmp,
767 
768  labelListList(0), // featurePointNormals,
769  labelListList(0), // featurePointEdges,
770  labelList(0) // regionEdges
771  );
772 
773  feMesh.write();
774 
775  feMesh.writeObj(feMesh.path()/sFeatFileName);
776 
777  {
778  // Write a featureEdgeMesh for backwards compatibility
779  featureEdgeMesh bfeMesh
780  (
781  IOobject
782  (
783  sFeatFileName + ".eMesh", // name
784  runTime.constant(), // instance
785  "triSurface",
786  runTime, // registry
789  false
790  ),
791  feMesh.points(),
792  feMesh.edges()
793  );
794 
795  Info<< nl << "Writing featureEdgeMesh to "
796  << bfeMesh.objectPath() << endl;
797 
798  bfeMesh.regIOobject::write();
799  }
800 
801  Info << "End\n" << endl;
802 
803  return 0;
804 }
805 
806 
807 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A class for handling file names.
Definition: fileName.H:69
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const double e
Elementary charge.
Definition: doubleFloat.H:78
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
void writeStats(Ostream &) const
Write some statistics.
Definition: triSurface.C:1067
Basic surface-surface intersection description. Constructed from two surfaces it creates a descriptio...
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const Type & second() const
Return second.
Definition: Pair.H:99
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const Field< PointType > & points() const
Return reference to global points.
static void noParallel()
Remove the parallel options.
Definition: argList.C:146
const Type & first() const
Return first.
Definition: Pair.H:87
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:154
void cross(FieldField< Field1, typename crossProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
virtual void movePoints(const pointField &)
Move points.
Definition: triSurface.C:788
Helper class to search on triSurface.
static scalarField minEdgeLength(const triSurface &surf)
Calculate min edge length for every surface point.
const stringList & args() const
Return arguments.
Definition: argListI.H:66
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
const pointField & points
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
A class for handling words, derived from string.
Definition: word.H:59
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
An STL-conforming hash table.
Definition: HashTable.H:61
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static scalar cosNormalAngleTol_
Angular closeness tolerance for treating normals as the same.
label start() const
Return start vertex label.
Definition: edgeI.H:81
static const char nl
Definition: Ostream.H:262
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.
Points connected by edges.
Definition: edgeMesh.H:69
A normal distribution model.
const Field< PointType > & faceNormals() const
Return face normals for patch.
vector vec(const pointField &) const
Return the vector (end - start)
Definition: edgeI.H:157
virtual bool write() const
Write using setting from DB.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:83
Triangulated surface description with patch information.
Definition: triSurface.H:65
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Holder of intersections of edges of a surface with another surface. Optionally shuffles around points...
Namespace for OpenFOAM.