surfaceIntersection.H
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 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 Class
25  Foam::surfaceIntersection
26 
27 Description
28  Basic surface-surface intersection description. Constructed from two
29  surfaces it creates a description of the intersection.
30 
31  The intersection information consists of the intersection line(s)
32  with new points, new edges between points (note that these edges and
33  points are on both surfaces) and various addressing from original
34  surface faces/edges to intersection and vice versa.
35 
36  Gets either precalculated intersection information or calculates it
37  itself.
38  Algorithm works by intersecting all edges of one surface with the other
39  surface and storing a reference from both faces (one on surface1, one on
40  surface 2) to the vertex. If the reference re-occurs we have the second
41  hit of both faces and an edge is created between the retrieved vertex and
42  the new one.
43 
44  Note: when doing intersecting itself uses intersection::planarTol() as a
45  fraction of
46  current edge length to determine if intersection is a point-touching one
47  instead of an edge-piercing action.
48 
49 SourceFiles
50  surfaceIntersection.C
51  surfaceIntersectionFuncs.C
52  surfaceIntersectionTemplates.C
53 
54 \*---------------------------------------------------------------------------*/
55 
56 #ifndef surfaceIntersection_H
57 #define surfaceIntersection_H
58 
59 #include "DynamicList.H"
60 #include "point.H"
61 #include "edge.H"
62 #include "labelPairLookup.H"
63 #include "typeInfo.H"
64 #include "edgeList.H"
65 #include "pointIndexHit.H"
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 namespace Foam
70 {
71 
72 // Forward declaration of classes
73 class triSurfaceSearch;
74 class triSurface;
75 class edgeIntersections;
76 
77 /*---------------------------------------------------------------------------*\
78  Class surfaceIntersection Declaration
79 \*---------------------------------------------------------------------------*/
80 
82 {
83  // Private data
84 
85  //- Newly introduced points.
87 
88  //- Newly introduced edges (are on both surfaces). Reference into
89  // cutPoints.
91 
92  //- From face on surf1 and face on surf2 to intersection point
93  // (label in cutPoints)
95 
96  //- From face on surf1 and face on surf2 to intersection edge
97  // (label in cutEdges)
99 
100  //- Edges on surf1 that are cut. From edge on surf1 to label in cutPoint
101  // If multiple cuts:sorted from edge.start to edge.end
103 
104  //- Edges on surf2 that are cut. From edge on surf2 to label in cutPoint
105  // If multiple cuts:sorted from edge.start to edge.end
107 
108 
109  // Private Member Functions
110 
111  //- Write point in obj format.
112  static void writeOBJ(const point& pt, Ostream& os);
113 
114  //- Write points and edges in obj format
115  static void writeOBJ
116  (
117  const List<point>&,
118  const List<edge>&,
119  Ostream&
120  );
121 
122  //- Transfer contents of List<DynamicList<..> > to List<List<..>>
123  template<class T>
124  static void transfer(List<DynamicList<T> >&, List<List<T> >&);
125 
126  //- Get minimum length of all edges connected to point
127  static scalar minEdgeLen(const triSurface& surf, const label pointI);
128 
129  //- Get edge label of edge between face vertices fp and fp+1
130  static label getEdge
131  (
132  const triSurface& surf,
133  const label faceI,
134  const label fp
135  );
136 
137  //- Remove duplicates from ordered dynamic list. Returns map from old
138  // to new (-1 if element removed)
139  static void removeDuplicates(const labelList& map, labelList& labels);
140 
141  //- Apply map to elements of a labelList
142  static void inlineRemap(const labelList& map, labelList& elems);
143 
144  // Remove all duplicate and degenerate elements. Return unique elements
145  // and map from old to new.
146  static edgeList filterEdges(const edgeList&, labelList& map);
147 
148  //- Remove all duplicate elements.
149  static labelList filterLabels(const labelList& elems, labelList& map);
150 
151  //- Do some checks if edge and face (resulting from hit)
152  // should not be considered. Returns true if can be discarded.
153  static bool excludeEdgeHit
154  (
155  const triSurface& surf,
156  const label edgeI,
157  const label faceI,
158  const scalar tol
159  );
160 
164  //static pointIndexHit faceEdgeIntersection
165  //(
166  // const triSurface&,
167  // const label hitFaceI,
168  //
169  // const vector& n,
170  // const point& eStart,
171  // const point& eEnd
172  //);
173 
174 
175  //- Debugging: Dump intersected edges to stream
177  (
178  const triSurface& surf,
179  const labelListList& edgeCutVerts,
180  Ostream& os
181  ) const;
182 
183  //- Detect if point close to edge of end. Returns -1: not close.
184  // 0:close (within startTol) to start, 1:close (within endTol) to end
185  static label classify
186  (
187  const scalar startTol,
188  const scalar endTol,
189  const point& p,
190  const edge& e,
191  const pointField& points
192  );
193 
194  //- Update reference between faceA and faceB. Updates facePairToVertex_
195  // (first occurrence of face pair) and facePairToEdge_ (second occ.)
196  void storeIntersection
197  (
198  const bool isFirstSurf,
199  const labelList& facesA,
200  const label faceB,
203  );
204 
205  //- Investigate pHit to whether is case of point hits point,
206  // point hits edge, point hits face or edge hits face.
207  void classifyHit
208  (
209  const triSurface& surf1,
210  const scalarField& surf1PointTol,
211  const triSurface& surf2,
212  const bool isFirstSurf,
213  const label edgeI,
214  const scalar tolDim,
215  const pointIndexHit& pHit,
216 
217  DynamicList<edge>& allCutEdges,
218  DynamicList<point>& allCutPoints,
219  List<DynamicList<label> >& surfEdgeCuts
220  );
221 
222  //- Cut edges of surf1 with surface 2.
223  void doCutEdges
224  (
225  const triSurface& surf1,
226  const triSurfaceSearch& querySurf2,
227  const bool isFirstSurf,
228  const bool isSelfIntersection,
229 
230  DynamicList<edge>& allCutEdges,
231  DynamicList<point>& allCutPoints,
232  List<DynamicList<label> >& surfEdgeCuts
233  );
234 
235 
236 public:
237 
238  ClassName("surfaceIntersection");
239 
240 
241  // Constructors
242 
243  //- Construct null
245 
246  //- Construct from precalculated intersection information.
247  // Advantage: intersection information is guaranteed to have no
248  // degenerate cuts.
250  (
251  const triSurface& surf1,
252  const edgeIntersections& intersections1,
253  const triSurface& surf2,
254  const edgeIntersections& intersections2
255  );
256 
257  //- Construct from two surfaces. Does all its own cutting.
258  // Has problems with degenerate cuts
260  (
261  const triSurfaceSearch& querySurf1,
262  const triSurfaceSearch& querySurf2
263  );
264 
265  //- Special: intersect surface with itself. Used to check for
266  // self-intersection.
267  surfaceIntersection(const triSurfaceSearch& querySurf1);
268 
269 
270  // Member Functions
271 
272  const pointField& cutPoints() const;
273 
274  const edgeList& cutEdges() const;
275 
276  const labelPairLookup& facePairToVertex() const;
277 
278  const labelPairLookup& facePairToEdge() const;
279 
280  //- Access either surf1EdgeCuts (isFirstSurface = true) or
281  // surf2EdgeCuts
282  const labelListList& edgeCuts(const bool) const;
283 
284  const labelListList& surf1EdgeCuts() const;
285 
286  const labelListList& surf2EdgeCuts() const;
287 };
288 
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 } // End namespace Foam
293 
294 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 
296 #ifdef NoRepository
298 #endif
299 
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 
302 #endif
303 
304 // ************************************************************************* //
Foam::surfaceIntersection::surf2EdgeCuts
const labelListList & surf2EdgeCuts() const
Definition: surfaceIntersection.C:1178
Foam::surfaceIntersection::facePairToVertex
const labelPairLookup & facePairToVertex() const
Definition: surfaceIntersection.C:1144
labelPairLookup.H
pointIndexHit.H
Foam::surfaceIntersection::facePairToEdge
const labelPairLookup & facePairToEdge() const
Definition: surfaceIntersection.C:1150
p
p
Definition: pEqn.H:62
typeInfo.H
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::surfaceIntersection::filterLabels
static labelList filterLabels(const labelList &elems, labelList &map)
Remove all duplicate elements.
Definition: surfaceIntersectionFuncs.C:223
point.H
Foam::surfaceIntersection::edgeCuts
const labelListList & edgeCuts(const bool) const
Access either surf1EdgeCuts (isFirstSurface = true) or.
Definition: surfaceIntersection.C:1157
Foam::surfaceIntersection::classify
static label classify(const scalar startTol, const scalar endTol, const point &p, const edge &e, const pointField &points)
Detect if point close to edge of end. Returns -1: not close.
Definition: surfaceIntersectionFuncs.C:304
Foam::surfaceIntersection::removeDuplicates
static void removeDuplicates(const labelList &map, labelList &labels)
Remove duplicates from ordered dynamic list. Returns map from old.
Definition: surfaceIntersectionFuncs.C:118
Foam::surfaceIntersection::cutPoints_
pointField cutPoints_
Newly introduced points.
Definition: surfaceIntersection.H:85
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::surfaceIntersection::doCutEdges
void doCutEdges(const triSurface &surf1, const triSurfaceSearch &querySurf2, const bool isFirstSurf, const bool isSelfIntersection, DynamicList< edge > &allCutEdges, DynamicList< point > &allCutPoints, List< DynamicList< label > > &surfEdgeCuts)
Cut edges of surf1 with surface 2.
Definition: surfaceIntersection.C:558
Foam::surfaceIntersection::filterEdges
static edgeList filterEdges(const edgeList &, labelList &map)
Definition: surfaceIntersectionFuncs.C:182
Foam::surfaceIntersection::surfaceIntersection
surfaceIntersection()
Construct null.
Definition: surfaceIntersection.C:699
Foam::surfaceIntersection::writeOBJ
static void writeOBJ(const point &pt, Ostream &os)
Write point in obj format.
Definition: surfaceIntersectionFuncs.C:37
Foam::surfaceIntersection::writeIntersectedEdges
void writeIntersectedEdges(const triSurface &surf, const labelListList &edgeCutVerts, Ostream &os) const
Debugging: Dump intersected edges to stream.
Definition: surfaceIntersectionFuncs.C:259
Foam::surfaceIntersection::ClassName
ClassName("surfaceIntersection")
Foam::triSurfaceSearch
Helper class to search on triSurface.
Definition: triSurfaceSearch.H:55
Foam::surfaceIntersection::getEdge
static label getEdge(const triSurface &surf, const label faceI, const label fp)
Get edge label of edge between face vertices fp and fp+1.
Definition: surfaceIntersectionFuncs.C:87
surfaceIntersectionTemplates.C
Foam::surfaceIntersection::surf1EdgeCuts
const labelListList & surf1EdgeCuts() const
Definition: surfaceIntersection.C:1172
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Foam::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::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::surfaceIntersection::surf2EdgeCuts_
labelListList surf2EdgeCuts_
Edges on surf2 that are cut. From edge on surf2 to label in cutPoint.
Definition: surfaceIntersection.H:105
Foam::surfaceIntersection::inlineRemap
static void inlineRemap(const labelList &map, labelList &elems)
Apply map to elements of a labelList.
Definition: surfaceIntersectionFuncs.C:167
Foam::surfaceIntersection::classifyHit
void classifyHit(const triSurface &surf1, const scalarField &surf1PointTol, const triSurface &surf2, const bool isFirstSurf, const label edgeI, const scalar tolDim, const pointIndexHit &pHit, DynamicList< edge > &allCutEdges, DynamicList< point > &allCutPoints, List< DynamicList< label > > &surfEdgeCuts)
Investigate pHit to whether is case of point hits point,.
Definition: surfaceIntersection.C:274
Foam::surfaceIntersection::cutEdges_
edgeList cutEdges_
Newly introduced edges (are on both surfaces). Reference into.
Definition: surfaceIntersection.H:89
Foam::surfaceIntersection::cutPoints
const pointField & cutPoints() const
Definition: surfaceIntersection.C:1132
edgeList.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::surfaceIntersection::cutEdges
const edgeList & cutEdges() const
Definition: surfaceIntersection.C:1138
edge.H
Foam::HashTable< label, FixedList< label, 2 >, FixedList< label, 2 >::Hash<> >
Foam::surfaceIntersection::transfer
static void transfer(List< DynamicList< T > > &, List< List< T > > &)
Transfer contents of List<DynamicList<..> > to List<List<..>>
Definition: surfaceIntersectionTemplates.C:33
Foam::surfaceIntersection::storeIntersection
void storeIntersection(const bool isFirstSurf, const labelList &facesA, const label faceB, DynamicList< edge > &, DynamicList< point > &)
Update reference between faceA and faceB. Updates facePairToVertex_.
Definition: surfaceIntersection.C:197
Foam::Vector< scalar >
Foam::surfaceIntersection::facePairToEdge_
labelPairLookup facePairToEdge_
From face on surf1 and face on surf2 to intersection edge.
Definition: surfaceIntersection.H:97
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
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::surfaceIntersection::surf1EdgeCuts_
labelListList surf1EdgeCuts_
Edges on surf1 that are cut. From edge on surf1 to label in cutPoint.
Definition: surfaceIntersection.H:101
DynamicList.H
Foam::surfaceIntersection::minEdgeLen
static scalar minEdgeLen(const triSurface &surf, const label pointI)
Get minimum length of all edges connected to point.
Definition: surfaceIntersectionFuncs.C:65
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::edgeIntersections
Holder of intersections of edges of a surface with another surface. Optionally shuffles around points...
Definition: edgeIntersections.H:60
Foam::surfaceIntersection::facePairToVertex_
labelPairLookup facePairToVertex_
From face on surf1 and face on surf2 to intersection point.
Definition: surfaceIntersection.H:93
Foam::surfaceIntersection::excludeEdgeHit
static bool excludeEdgeHit(const triSurface &surf, const label edgeI, const label faceI, const scalar tol)
Do some checks if edge and face (resulting from hit)
Definition: surfaceIntersection.C:52