surfaceIntersectionFuncs.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "surfaceIntersection.H"
30 #include "triSurface.H"
31 #include "triSurfaceSearch.H"
32 #include "edgeHashes.H"
33 #include "labelPairHashes.H"
34 #include "OFstream.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // file-scope
42 // Write points in obj format
43 static void writeObjPoints(const UList<point>& pts, Ostream& os)
44 {
45  for (const point& pt : pts)
46  {
47  os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
48  }
49 }
50 
51 } // End namespace Foam
52 
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 
57 // Get minimum length of all edges connected to point
58 Foam::scalar Foam::surfaceIntersection::minEdgeLen
59 (
60  const triSurface& surf,
61  const label pointi
62 )
63 {
64  const labelList& pEdges = surf.pointEdges()[pointi];
65 
66  scalar minLen = GREAT;
67 
68  forAll(pEdges, pEdgeI)
69  {
70  const edge& e = surf.edges()[pEdges[pEdgeI]];
71 
72  minLen = min(minLen, e.mag(surf.localPoints()));
73  }
74 
75  return minLen;
76 }
77 
78 
79 // Get edge between fp and fp+1 on facei.
80 Foam::label Foam::surfaceIntersection::getEdge
81 (
82  const triSurface& surf,
83  const label facei,
84  const label fp
85 )
86 {
87  const edge faceEdge = surf.localFaces()[facei].faceEdge(fp);
88 
89  const labelList& eLabels = surf.faceEdges()[facei];
90 
91  forAll(eLabels, eI)
92  {
93  const label edgeI = eLabels[eI];
94 
95  if (surf.edges()[edgeI] == faceEdge)
96  {
97  return edgeI;
98  }
99  }
100 
102  << "Problem:: Cannot find edge with vertices " << faceEdge
103  << " in face " << facei
104  << abort(FatalError);
105 
106  return -1;
107 }
108 
109 
110 // Given a map remove all consecutive duplicate elements.
111 void Foam::surfaceIntersection::removeDuplicates
112 (
113  const labelList& map,
114  labelList& elems
115 )
116 {
117  bool hasDuplicate = false;
118 
119  label prevVertI = -1;
120 
121  forAll(elems, elemI)
122  {
123  label newVertI = map[elems[elemI]];
124 
125  if (newVertI == prevVertI)
126  {
127  hasDuplicate = true;
128 
129  break;
130  }
131  prevVertI = newVertI;
132  }
133 
134  if (hasDuplicate)
135  {
136  // Create copy
137  labelList oldElems(elems);
138 
139  label elemI = 0;
140 
141  // Insert first
142  elems[elemI++] = map[oldElems[0]];
143 
144  for (label vertI = 1; vertI < oldElems.size(); vertI++)
145  {
146  // Insert others only if they differ from one before
147  label newVertI = map[oldElems[vertI]];
148 
149  if (newVertI != elems.last())
150  {
151  elems[elemI++] = newVertI;
152  }
153  }
154  elems.setSize(elemI);
155  }
156 }
157 
158 
159 // Remove all duplicate and degenerate elements. Return unique elements and
160 // map from old to new.
161 Foam::edgeList Foam::surfaceIntersection::filterEdges
162 (
163  const edgeList& edges,
164  labelList& map
165 )
166 {
167  edgeHashSet uniqueEdges(10*edges.size());
168 
169  edgeList newEdges(edges.size());
170 
171  map.setSize(edges.size());
172  map = -1;
173 
174  label newEdgeI = 0;
175 
176  forAll(edges, edgeI)
177  {
178  const edge& e = edges[edgeI];
179 
180  if ((e.start() != e.end()) && uniqueEdges.insert(e))
181  {
182  // Edge is non-degenerate and not yet seen.
183 
184  map[edgeI] = newEdgeI;
185 
186  newEdges[newEdgeI++] = e;
187  }
188  }
189 
190  newEdges.setSize(newEdgeI);
191 
192  return newEdges;
193 }
194 
195 
196 // Remove all duplicate elements.
197 Foam::labelList Foam::surfaceIntersection::filterLabels
198 (
199  const labelList& elems,
200  labelList& map
201 )
202 {
203  labelHashSet uniqueElems(10*elems.size());
204 
205  labelList newElems(elems.size());
206 
207  map.setSize(elems.size());
208  map = -1;
209 
210  label newElemI = 0;
211 
212  forAll(elems, elemI)
213  {
214  label elem = elems[elemI];
215 
216  if (uniqueElems.insert(elem))
217  {
218  // First time elem is seen
219 
220  map[elemI] = newElemI;
221 
222  newElems[newElemI++] = elem;
223  }
224  }
225 
226  newElems.setSize(newElemI);
227 
228  return newElems;
229 }
230 
231 
232 void Foam::surfaceIntersection::writeIntersectedEdges
233 (
234  const triSurface& surf,
235  const labelListList& edgeCutVerts,
236  Ostream& os
237 ) const
238 {
239  // Dump all points (surface followed by cutPoints)
240  const pointField& pts = surf.localPoints();
241 
242  writeObjPoints(pts, os);
243  writeObjPoints(cutPoints(), os);
244 
245  forAll(edgeCutVerts, edgeI)
246  {
247  const labelList& extraVerts = edgeCutVerts[edgeI];
248 
249  if (extraVerts.size())
250  {
251  const edge& e = surf.edges()[edgeI];
252 
253  // Start of original edge to first extra point
254  os << "l " << e.start()+1 << ' '
255  << extraVerts[0] + surf.nPoints() + 1 << nl;
256 
257  for (label i = 1; i < extraVerts.size(); i++)
258  {
259  os << "l " << extraVerts[i-1] + surf.nPoints() + 1 << ' '
260  << extraVerts[i] + surf.nPoints() + 1 << nl;
261  }
262 
263  os << "l " << extraVerts.last() + surf.nPoints() + 1
264  << ' ' << e.end()+1 << nl;
265  }
266  }
267 }
268 
269 
270 // Return 0 (p close to start), 1(close to end) or -1.
271 Foam::label Foam::surfaceIntersection::classify
272 (
273  const scalar startTol,
274  const scalar endTol,
275  const point& p,
276  const edge& e,
277  const UList<point>& points
278 )
279 {
280  if (mag(p - points[e.start()]) < startTol)
281  {
282  return 0;
283  }
284  else if (mag(p - points[e.end()]) < endTol)
285  {
286  return 1;
287  }
288 
289  return -1;
290 }
291 
292 
293 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::edgeList
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
triSurface.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
surfaceIntersection.H
labelPairHashes.H
A HashTable to objects of type <T> with a labelPair key. The hashing is based on labelPair (FixedList...
Foam::FatalError
error FatalError
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::writeObjPoints
static void writeObjPoints(const UList< point > &pts, Ostream &os)
Definition: edgeSurface.C:42
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
edgeHashes.H
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::Vector< scalar >
Foam::List< edge >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
triSurfaceSearch.H
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
Definition: HashSet.H:409
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::edgeHashSet
HashSet< edge, Hash< edge > > edgeHashSet
A HashSet with edge for its key.
Definition: edgeHashes.H:51