1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
26 Foam::surfaceIntersection
29 Basic surface-surface intersection description. Constructed from two
30 surfaces it creates a description of the intersection.
32 The intersection information consists of the intersection line(s)
33 with new points, new edges between points (note that these edges and
34 points are on both surfaces) and various addressing from original
35 surface faces/edges to intersection and vice versa.
37 Gets either precalculated intersection information or calculates it
39 Algorithm works by intersecting all edges of one surface with the other
40 surface and storing a reference from both faces (one on surface1, one on
41 surface 2) to the vertex. If the reference re-occurs we have the second
42 hit of both faces and an edge is created between the retrieved vertex and
45 Note: when doing intersecting itself uses intersection::planarTol() as a
47 current edge length to determine if intersection is a point-touching one
48 instead of an edge-piercing action.
52 surfaceIntersectionFuncs.C
53 surfaceIntersectionTemplates.C
55 \*---------------------------------------------------------------------------*/
57 #ifndef surfaceIntersection_H
58 #define surfaceIntersection_H
60 #include "DynamicList.H"
63 #include "labelPairLookup.H"
66 #include "pointIndexHit.H"
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 // Forward declaration of classes
74 class triSurfaceSearch;
76 class edgeIntersections;
78 /*---------------------------------------------------------------------------*\
79 Class surfaceIntersection Declaration
80 \*---------------------------------------------------------------------------*/
82 class surfaceIntersection
86 //- Newly introduced points.
87 pointField cutPoints_;
89 //- Newly introduced edges (are on both surfaces). Reference into
93 //- From face on surf1 and face on surf2 to intersection point
94 // (label in cutPoints)
95 labelPairLookup facePairToVertex_;
97 //- From face on surf1 and face on surf2 to intersection edge
98 // (label in cutEdges)
99 labelPairLookup facePairToEdge_;
101 //- Edges on surf1 that are cut. From edge on surf1 to label in cutPoint
102 // If multiple cuts:sorted from edge.start to edge.end
103 labelListList surf1EdgeCuts_;
105 //- Edges on surf2 that are cut. From edge on surf2 to label in cutPoint
106 // If multiple cuts:sorted from edge.start to edge.end
107 labelListList surf2EdgeCuts_;
110 // Private Member Functions
112 //- Write point in obj format.
113 static void writeOBJ(const point& pt, Ostream& os);
115 //- Write points and edges in obj format
123 //- Transfer contents of List<DynamicList<..> > to List<List<..>>
125 static void transfer(List<DynamicList<T> >&, List<List<T> >&);
127 //- Get minimum length of all edges connected to point
128 static scalar minEdgeLen(const triSurface& surf, const label pointI);
130 //- Get edge label of edge between face vertices fp and fp+1
133 const triSurface& surf,
138 //- Remove duplicates from ordered dynamic list. Returns map from old
139 // to new (-1 if element removed)
140 static void removeDuplicates(const labelList& map, labelList& labels);
142 //- Apply map to elements of a labelList
143 static void inlineRemap(const labelList& map, labelList& elems);
145 // Remove all duplicate and degenerate elements. Return unique elements
146 // and map from old to new.
147 static edgeList filterEdges(const edgeList&, labelList& map);
149 //- Remove all duplicate elements.
150 static labelList filterLabels(const labelList& elems, labelList& map);
152 //- Do some checks if edge and face (resulting from hit)
153 // should not be considered. Returns true if can be discarded.
154 static bool excludeEdgeHit
156 const triSurface& surf,
162 ////- Given edge (eStart - eEnd) and normal direction construct plane
163 //// and intersect all edges of hitFace with it.
164 //// Return the edge and coordinate of hit.
165 //static pointIndexHit faceEdgeIntersection
167 // const triSurface&,
168 // const label hitFaceI,
171 // const point& eStart,
176 //- Debugging: Dump intersected edges to stream
177 void writeIntersectedEdges
179 const triSurface& surf,
180 const labelListList& edgeCutVerts,
184 //- Detect if point close to edge of end. Returns -1: not close.
185 // 0:close (within startTol) to start, 1:close (within endTol) to end
186 static label classify
188 const scalar startTol,
192 const pointField& points
195 //- Update reference between faceA and faceB. Updates facePairToVertex_
196 // (first occurrence of face pair) and facePairToEdge_ (second occ.)
197 void storeIntersection
199 const bool isFirstSurf,
200 const labelList& facesA,
206 //- Investigate pHit to whether is case of point hits point,
207 // point hits edge, point hits face or edge hits face.
210 const triSurface& surf1,
211 const scalarField& surf1PointTol,
212 const triSurface& surf2,
213 const bool isFirstSurf,
216 const pointIndexHit& pHit,
218 DynamicList<edge>& allCutEdges,
219 DynamicList<point>& allCutPoints,
220 List<DynamicList<label> >& surfEdgeCuts
223 //- Cut edges of surf1 with surface 2.
226 const triSurface& surf1,
227 const triSurfaceSearch& querySurf2,
228 const bool isFirstSurf,
229 const bool isSelfIntersection,
231 DynamicList<edge>& allCutEdges,
232 DynamicList<point>& allCutPoints,
233 List<DynamicList<label> >& surfEdgeCuts
238 ClassName("surfaceIntersection");
244 surfaceIntersection();
246 //- Construct from precalculated intersection information.
247 // Advantage: intersection information is guaranteed to have no
251 const triSurface& surf1,
252 const edgeIntersections& intersections1,
253 const triSurface& surf2,
254 const edgeIntersections& intersections2
257 //- Construct from two surfaces. Does all its own cutting.
258 // Has problems with degenerate cuts
261 const triSurfaceSearch& querySurf1,
262 const triSurfaceSearch& querySurf2
265 //- Special: intersect surface with itself. Used to check for
266 // self-intersection.
267 surfaceIntersection(const triSurfaceSearch& querySurf1);
272 const pointField& cutPoints() const;
274 const edgeList& cutEdges() const;
276 //const labelPairLookup& facePairToVertex() const;
278 const labelPairLookup& facePairToEdge() const;
280 //- Access either surf1EdgeCuts (isFirstSurface = true) or
282 const labelListList& edgeCuts(const bool) const;
284 const labelListList& surf1EdgeCuts() const;
286 const labelListList& surf2EdgeCuts() const;
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293 } // End namespace Foam
295 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
298 # include "surfaceIntersectionTemplates.C"
301 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305 // ************************************************************************* //