orientedSurface.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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "orientedSurface.H"
29 #include "triSurfaceTools.H"
30 #include "triSurfaceSearch.H"
31 #include "treeBoundBox.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 defineTypeNameAndDebug(orientedSurface, 0);
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 // Return true if edge is used in opposite order in faces
44 bool Foam::orientedSurface::consistentEdge
45 (
46  const edge& e,
47  const triSurface::FaceType& f0,
48  const triSurface::FaceType& f1
49 )
50 {
51  return (f0.edgeDirection(e) > 0) ^ (f1.edgeDirection(e) > 0);
52 }
53 
54 
55 Foam::labelList Foam::orientedSurface::faceToEdge
56 (
57  const triSurface& s,
58  const labelList& changedFaces
59 )
60 {
61  labelList changedEdges(3*changedFaces.size());
62  label changedI = 0;
63 
64  forAll(changedFaces, i)
65  {
66  const labelList& fEdges = s.faceEdges()[changedFaces[i]];
67 
68  forAll(fEdges, j)
69  {
70  changedEdges[changedI++] = fEdges[j];
71  }
72  }
73  changedEdges.setSize(changedI);
74 
75  return changedEdges;
76 }
77 
78 
79 Foam::labelList Foam::orientedSurface::edgeToFace
80 (
81  const triSurface& s,
82  const labelList& changedEdges,
83  labelList& flip
84 )
85 {
86  labelList changedFaces(2*changedEdges.size());
87  label changedI = 0;
88 
89  forAll(changedEdges, i)
90  {
91  label edgeI = changedEdges[i];
92 
93  const labelList& eFaces = s.edgeFaces()[edgeI];
94 
95  if (eFaces.size() < 2)
96  {
97  // Do nothing, faces was already visited.
98  }
99  else if (eFaces.size() == 2)
100  {
101  label face0 = eFaces[0];
102  label face1 = eFaces[1];
103 
104  const triSurface::FaceType& f0 = s.localFaces()[face0];
105  const triSurface::FaceType& f1 = s.localFaces()[face1];
106 
107  if (flip[face0] == UNVISITED)
108  {
109  if (flip[face1] == UNVISITED)
110  {
112  << abort(FatalError);
113  }
114  else
115  {
116  // Face1 has a flip state, face0 hasn't
117  if (consistentEdge(s.edges()[edgeI], f0, f1))
118  {
119  // Take over flip status
120  flip[face0] = (flip[face1] == FLIP ? FLIP : NOFLIP);
121  }
122  else
123  {
124  // Invert
125  flip[face0] = (flip[face1] == FLIP ? NOFLIP : FLIP);
126  }
127  changedFaces[changedI++] = face0;
128  }
129  }
130  else
131  {
132  if (flip[face1] == UNVISITED)
133  {
134  // Face0 has a flip state, face1 hasn't
135  if (consistentEdge(s.edges()[edgeI], f0, f1))
136  {
137  flip[face1] = (flip[face0] == FLIP ? FLIP : NOFLIP);
138  }
139  else
140  {
141  flip[face1] = (flip[face0] == FLIP ? NOFLIP : FLIP);
142  }
143  changedFaces[changedI++] = face1;
144  }
145  }
146  }
147  else
148  {
149  // Multiply connected. Do what?
150  }
151  }
152  changedFaces.setSize(changedI);
153 
154  return changedFaces;
155 }
156 
157 
158 void Foam::orientedSurface::walkSurface
159 (
160  const triSurface& s,
161  const label startFacei,
162  labelList& flipState
163 )
164 {
165  // List of faces that were changed in the last iteration.
166  labelList changedFaces(1, startFacei);
167  // List of edges that were changed in the last iteration.
168  labelList changedEdges;
169 
170  while (true)
171  {
172  changedEdges = faceToEdge(s, changedFaces);
173 
174  if (changedEdges.empty())
175  {
176  break;
177  }
178 
179  changedFaces = edgeToFace(s, changedEdges, flipState);
180 
181  if (changedFaces.empty())
182  {
183  break;
184  }
185  }
186 }
187 
188 
189 void Foam::orientedSurface::propagateOrientation
190 (
191  const triSurface& s,
192  const point& samplePoint,
193  const bool orientOutside,
194  const label nearestFacei,
195  const point& nearestPt,
196  labelList& flipState
197 )
198 {
199  //
200  // Determine orientation to normal on nearest face
201  //
203  (
204  s,
205  samplePoint,
206  nearestFacei
207  );
208 
209  if (side == triSurfaceTools::UNKNOWN)
210  {
211  // Non-closed surface. Do what? For now behave as if no flipping
212  // necessary
213  flipState[nearestFacei] = NOFLIP;
214  }
215  else if ((side == triSurfaceTools::OUTSIDE) == orientOutside)
216  {
217  // outside & orientOutside or inside & !orientOutside
218  // Normals on surface pointing correctly. No need to flip normals
219  flipState[nearestFacei] = NOFLIP;
220  }
221  else
222  {
223  // Need to flip normals.
224  flipState[nearestFacei] = FLIP;
225  }
226 
227  if (debug)
228  {
229  vector n = triSurfaceTools::surfaceNormal(s, nearestFacei, nearestPt);
230 
231  Pout<< "orientedSurface::propagateOrientation : starting face"
232  << " orientation:" << nl
233  << " for samplePoint:" << samplePoint << nl
234  << " starting from point:" << nearestPt << nl
235  << " on face:" << nearestFacei << nl
236  << " with normal:" << n << nl
237  << " decided side:" << label(side)
238  << endl;
239  }
240 
241  // Walk the surface from nearestFacei, changing the flipstate.
242  walkSurface(s, nearestFacei, flipState);
243 }
244 
245 
246 // Find side for zoneI only by counting the number of intersections. Determines
247 // if face is oriented consistent with outwards pointing normals.
248 void Foam::orientedSurface::findZoneSide
249 (
250  const triSurfaceSearch& surfSearches,
251  const labelList& faceZone,
252  const label zoneI,
253  const point& outsidePoint,
254  label& zoneFacei,
255  bool& isOutside
256 )
257 {
258  const triSurface& s = surfSearches.surface();
259 
260  zoneFacei = -1;
261  isOutside = false;
262 
263  pointField start(1, outsidePoint);
264  List<List<pointIndexHit>> hits(1, List<pointIndexHit>());
265 
266  forAll(faceZone, facei)
267  {
268  if (faceZone[facei] == zoneI)
269  {
270  const point& fc = s.faceCentres()[facei];
271  const vector& n = s.faceNormals()[facei];
272 
273  const vector d = fc - outsidePoint;
274  const scalar magD = mag(d);
275 
276  // Check if normal different enough to decide upon
277  if (magD > SMALL && (mag(n & d/magD) > 1e-6))
278  {
279  pointField end(1, fc + d);
280 
281  //Info<< "Zone " << zoneI << " : Shooting ray"
282  // << " from " << outsidePoint
283  // << " to " << end
284  // << " to pierce triangle " << facei
285  // << " with centre " << fc << endl;
286 
287 
288  surfSearches.findLineAll(start, end, hits);
289 
290  label zoneIndex = -1;
291  forAll(hits[0], i)
292  {
293  if (hits[0][i].index() == facei)
294  {
295  zoneIndex = i;
296  break;
297  }
298  }
299 
300  if (zoneIndex != -1)
301  {
302  zoneFacei = facei;
303 
304  if ((zoneIndex%2) == 0)
305  {
306  // Odd number of intersections. Check if normal points
307  // in direction of ray
308  isOutside = ((n & d) < 0);
309  }
310  else
311  {
312  isOutside = ((n & d) > 0);
313  }
314  break;
315  }
316  }
317  }
318  }
319 }
320 
321 
322 bool Foam::orientedSurface::flipSurface
323 (
324  triSurface& s,
325  const labelList& flipState
326 )
327 {
328  bool hasFlipped = false;
329 
330  // Flip tris in s
331  forAll(flipState, facei)
332  {
333  if (flipState[facei] == UNVISITED)
334  {
336  << "unvisited face " << facei
337  << abort(FatalError);
338  }
339  else if (flipState[facei] == FLIP)
340  {
341  labelledTri& tri = s[facei];
342 
343  label tmp = tri[0];
344 
345  tri[0] = tri[1];
346  tri[1] = tmp;
347 
348  hasFlipped = true;
349  }
350  }
351  // Recalculate normals
352  if (hasFlipped)
353  {
354  s.clearOut();
355  }
356  return hasFlipped;
357 }
358 
359 
360 bool Foam::orientedSurface::orientConsistent(triSurface& s)
361 {
362  bool anyFlipped = false;
363 
364  // Do initial flipping to make triangles consistent. Otherwise if the
365  // nearest is e.g. on an edge inbetween inconsistent triangles it might
366  // make the wrong decision.
367  if (s.size() > 0)
368  {
369  // Whether face has to be flipped.
370  // UNVISITED: unvisited
371  // NOFLIP: no need to flip
372  // FLIP: need to flip
373  labelList flipState(s.size(), UNVISITED);
374 
375  label facei = 0;
376  while (true)
377  {
378  label startFacei = -1;
379  while (facei < s.size())
380  {
381  if (flipState[facei] == UNVISITED)
382  {
383  startFacei = facei;
384  break;
385  }
386  facei++;
387  }
388 
389  if (startFacei == -1)
390  {
391  break;
392  }
393 
394  flipState[startFacei] = NOFLIP;
395  walkSurface(s, startFacei, flipState);
396  }
397 
398  anyFlipped = flipSurface(s, flipState);
399  }
400  return anyFlipped;
401 }
402 
403 
404 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
405 
406 // Null constructor
408 :
409  triSurface()
410 {}
411 
412 
413 // Construct from surface and point which defines outside
415 (
416  const triSurface& surf,
417  const point& samplePoint,
418  const bool orientOutside
419 )
420 :
421  triSurface(surf)
422 {
423  orient(*this, samplePoint, orientOutside);
424 }
425 
426 
427 // Construct from surface. Calculate outside point.
429 (
430  const triSurface& surf,
431  const bool orientOutside
432 )
433 :
434  triSurface(surf)
435 {
436  // BoundBox calculation without localPoints
437  treeBoundBox bb(surf.points(), surf.meshPoints());
438 
439  point outsidePoint = bb.max() + bb.span();
440 
441  orient(*this, outsidePoint, orientOutside);
442 }
443 
444 
445 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
446 
448 (
449  triSurface& s,
450  const point& samplePoint,
451  const bool orientOutside
452 )
453 {
454  // Do initial flipping to make triangles consistent. Otherwise if the
455  // nearest is e.g. on an edge inbetween inconsistent triangles it might
456  // make the wrong decision.
457  bool topoFlipped = orientConsistent(s);
458 
459 
460  // Whether face has to be flipped.
461  // UNVISITED: unvisited
462  // NOFLIP: no need to flip
463  // FLIP: need to flip
464  labelList flipState(s.size(), UNVISITED);
465 
466 
467  while (true)
468  {
469  // Linear search for nearest unvisited point on surface.
470 
471  scalar minDist = GREAT;
472  point minPoint;
473  label minFacei = -1;
474 
475  forAll(s, facei)
476  {
477  if (flipState[facei] == UNVISITED)
478  {
479  pointHit curHit =
480  s[facei].nearestPoint(samplePoint, s.points());
481 
482  if (curHit.distance() < minDist)
483  {
484  minDist = curHit.distance();
485  minPoint = curHit.rawPoint();
486  minFacei = facei;
487  }
488  }
489  }
490 
491  // Did we find anything?
492  if (minFacei == -1)
493  {
494  break;
495  }
496 
497  // From this nearest face see if needs to be flipped and then
498  // go outwards.
499  propagateOrientation
500  (
501  s,
502  samplePoint,
503  orientOutside,
504  minFacei,
505  minPoint,
506  flipState
507  );
508  }
509 
510  // Now finally flip triangles according to flipState.
511  bool geomFlipped = flipSurface(s, flipState);
512 
513  return topoFlipped || geomFlipped;
514 }
515 
516 
518 (
519  triSurface& s,
520  const triSurfaceSearch& querySurf,
521  const point& samplePoint,
522  const bool orientOutside
523 )
524 {
525  // Do initial flipping to make triangles consistent. Otherwise if the
526  // nearest is e.g. on an edge inbetween inconsistent triangles it might
527  // make the wrong decision.
528  bool topoFlipped = orientConsistent(s);
529 
530  // Determine disconnected parts of surface
531  boolList borderEdge(s.nEdges(), false);
532  forAll(s.edgeFaces(), edgeI)
533  {
534  if (s.edgeFaces()[edgeI].size() != 2)
535  {
536  borderEdge[edgeI] = true;
537  }
538  }
540  label nZones = s.markZones(borderEdge, faceZone);
541 
542  // Check intersection with one face per zone.
543 
544  labelList flipState(s.size(), UNVISITED);
545  for (label zoneI = 0; zoneI < nZones; zoneI++)
546  {
547  label zoneFacei = -1;
548  bool isOutside;
549  findZoneSide
550  (
551  querySurf,
552  faceZone,
553  zoneI,
554  samplePoint,
555 
556  zoneFacei,
557  isOutside
558  );
559 
560  if (isOutside == orientOutside)
561  {
562  flipState[zoneFacei] = NOFLIP;
563  }
564  else
565  {
566  flipState[zoneFacei] = FLIP;
567  }
568  walkSurface(s, zoneFacei, flipState);
569  }
570 
571  // Now finally flip triangles according to flipState.
572  bool geomFlipped = flipSurface(s, flipState);
573 
574  return topoFlipped || geomFlipped;
575 }
576 
577 
578 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
nZones
label nZones
Definition: interpolatedFaces.H:24
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatch.H:300
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::PointHit< point >
Foam::triSurfaceTools::surfaceNormal
static vector surfaceNormal(const triSurface &surf, const label nearestFacei, const point &nearestPt)
Triangle (unit) normal. If nearest point to triangle on edge use.
Definition: triSurfaceTools.C:1955
Foam::orientedSurface::orientedSurface
orientedSurface()
Construct null.
Definition: orientedSurface.C:407
Foam::PointHit::rawPoint
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:162
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:87
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::triSurfaceSearch
Helper class to search on triSurface.
Definition: triSurfaceSearch.H:58
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::PointHit::distance
scalar distance() const
Return distance to hit.
Definition: PointHit.H:143
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:62
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:70
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
treeBoundBox.H
Foam::PrimitivePatch< labelledTri, ::Foam::List, pointField, point >::FaceType
labelledTri FaceType
Definition: PrimitivePatch.H:100
Foam::FatalError
error FatalError
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:115
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::triSurfaceTools::UNKNOWN
Definition: triSurfaceTools.H:452
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::triSurfaceTools::surfaceSide
static sideType surfaceSide(const triSurface &surf, const point &sample, const label nearestFacei)
Given nearest point (to sample) on surface determines which side.
Definition: triSurfaceTools.C:2036
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::VectorSpace::max
static const Form max
Definition: VectorSpace.H:117
Foam::Vector< scalar >
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
Foam::triSurfaceTools::OUTSIDE
Definition: triSurfaceTools.H:454
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
triSurfaceSearch.H
Foam::orientedSurface::orient
static bool orient(triSurface &, const point &, const bool orientOutside)
Flip faces such that normals are consistent with point:
Definition: orientedSurface.C:448
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:418
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::triSurfaceTools::sideType
sideType
On which side of surface.
Definition: triSurfaceTools.H:450
triSurfaceTools.H
orientedSurface.H