edgeIntersections.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 \*---------------------------------------------------------------------------*/
25 
26 #include "edgeIntersections.H"
27 #include "triSurfaceSearch.H"
28 #include "labelPairLookup.H"
29 #include "OFstream.H"
30 #include "HashSet.H"
31 #include "triSurface.H"
32 #include "pointIndexHit.H"
33 #include "treeDataTriSurface.H"
34 #include "indexedOctree.H"
35 #include "meshTools.H"
36 #include "plane.H"
37 #include "Random.H"
38 #include "unitConversion.H"
39 #include "treeBoundBox.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 defineTypeNameAndDebug(edgeIntersections, 0);
46 
47 scalar edgeIntersections::alignedCos_ = cos(degToRad(89.0));
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 void Foam::edgeIntersections::checkEdges(const triSurface& surf)
54 {
55  const pointField& localPoints = surf.localPoints();
56  const edgeList& edges = surf.edges();
57  const labelListList& edgeFaces = surf.edgeFaces();
58 
59  treeBoundBox bb(localPoints);
60 
61  scalar minSize = small * bb.minDim();
62 
63  forAll(edges, edgeI)
64  {
65  const edge& e = edges[edgeI];
66 
67  scalar eMag = e.mag(localPoints);
68 
69  if (eMag < minSize)
70  {
72  << "Edge " << edgeI << " vertices " << e
73  << " coords:" << localPoints[e[0]] << ' '
74  << localPoints[e[1]] << " is very small compared to bounding"
75  << " box dimensions " << bb << endl
76  << "This might lead to problems in intersection"
77  << endl;
78  }
79 
80  if (edgeFaces[edgeI].size() == 1)
81  {
83  << "Edge " << edgeI << " vertices " << e
84  << " coords:" << localPoints[e[0]] << ' '
85  << localPoints[e[1]] << " has only one face connected to it:"
86  << edgeFaces[edgeI] << endl
87  << "This might lead to problems in intersection"
88  << endl;
89  }
90  }
91 }
92 
93 
94 // Update intersections for selected edges.
95 void Foam::edgeIntersections::intersectEdges
96 (
97  const triSurface& surf1,
98  const pointField& points1, // surf1 meshPoints (not localPoints!)
99  const triSurfaceSearch& querySurf2,
100  const scalarField& surf1PointTol, // surf1 tolerance per point
101  const labelList& edgeLabels
102 )
103 {
104  const triSurface& surf2 = querySurf2.surface();
105  const vectorField& normals2 = surf2.faceNormals();
106 
107  const labelList& meshPoints = surf1.meshPoints();
108 
109  if (debug)
110  {
111  Pout<< "Calculating intersection of " << edgeLabels.size() << " edges"
112  << " out of " << surf1.nEdges() << " with " << surf2.size()
113  << " triangles ..." << endl;
114  }
115 
116  pointField start(edgeLabels.size());
117  pointField end(edgeLabels.size());
118  vectorField edgeDirs(edgeLabels.size());
119 
120  // Go through all edges, calculate intersections
121  forAll(edgeLabels, i)
122  {
123  label edgeI = edgeLabels[i];
124 
125  if (debug)// && (i % 1000 == 0))
126  {
127  Pout<< "Intersecting edge " << edgeI << " with surface" << endl;
128  }
129 
130  const edge& e = surf1.edges()[edgeI];
131 
132  const point& pStart = points1[meshPoints[e.start()]];
133  const point& pEnd = points1[meshPoints[e.end()]];
134 
135  const vector eVec(pEnd - pStart);
136  const vector n(eVec/(mag(eVec) + vSmall));
137 
138  // Start tracking somewhat before pStart and up to somewhat after p1.
139  // Note that tolerances here are smaller than those used to classify
140  // hit below.
141  // This will cause this hit to be marked as degenerate and resolved
142  // later on.
143  start[i] = pStart - 0.5*surf1PointTol[e[0]]*n;
144  end[i] = pEnd + 0.5*surf1PointTol[e[1]]*n;
145 
146  edgeDirs[i] = n;
147  }
148 
149  List<List<pointIndexHit>> edgeIntersections;
150  querySurf2.findLineAll
151  (
152  start,
153  end,
154  edgeIntersections
155  );
156 
157  label nHits = 0;
158 
159  // Classify the hits
160  forAll(edgeLabels, i)
161  {
162  const label edgeI = edgeLabels[i];
163 
164  labelList& intersectionTypes = classification_[edgeI];
165  intersectionTypes.setSize(edgeIntersections[i].size(), -1);
166 
167  this->operator[](edgeI).transfer(edgeIntersections[i]);
168 
169  forAll(intersectionTypes, hitI)
170  {
171  const pointIndexHit& pHit = this->operator[](edgeI)[hitI];
172  label& hitType = intersectionTypes[hitI];
173 
174  if (!pHit.hit())
175  {
176  continue;
177  }
178 
179  const edge& e = surf1.edges()[edgeI];
180 
181  // Classify point on surface1 edge.
182  if (mag(pHit.hitPoint() - start[i]) < surf1PointTol[e[0]])
183  {
184  // Intersection is close to edge start
185  hitType = 0;
186  }
187  else if (mag(pHit.hitPoint() - end[i]) < surf1PointTol[e[1]])
188  {
189  // Intersection is close to edge end
190  hitType = 1;
191  }
192  else if (mag(edgeDirs[i] & normals2[pHit.index()]) < alignedCos_)
193  {
194  // Edge is almost coplanar with the face
195 
196  Pout<< "Flat angle edge:" << edgeI
197  << " face:" << pHit.index()
198  << " cos:" << mag(edgeDirs[i] & normals2[pHit.index()])
199  << endl;
200 
201  hitType = 2;
202  }
203 
204  if (debug)
205  {
206  Info<< " hit " << pHit << " classify = " << hitType << endl;
207  }
208 
209  nHits++;
210  }
211  }
212 
213  if (debug)
214  {
215  Pout<< "Found " << nHits << " intersections of edges with surface ..."
216  << endl;
217  }
218 }
219 
220 
221 // If edgeI intersections are close to endpoint of edge shift endpoints
222 // slightly along edge
223 // Updates
224 // - points1 with new endpoint position
225 // - affectedEdges with all edges affected by moving the point
226 // Returns true if changed anything.
227 bool Foam::edgeIntersections::inlinePerturb
228 (
229  const triSurface& surf1,
230  const scalarField& surf1PointTol, // surf1 tolerance per point
231  const label edgeI,
232  Random& rndGen,
233  pointField& points1,
234  boolList& affectedEdges
235 ) const
236 {
237  bool hasPerturbed = false;
238 
239  // Check if edge close to endpoint. Note that we only have to check
240  // the intersection closest to the edge endpoints (i.e. first and last in
241  // edgeEdgs)
242 
243  const labelList& edgeEnds = classification_[edgeI];
244 
245  if (edgeEnds.size())
246  {
247  bool perturbStart = false;
248  bool perturbEnd = false;
249 
250  // Check first intersection.
251  if (edgeEnds.first() == 0)
252  {
253  perturbStart = true;
254  }
255 
256  if (edgeEnds.last() == 1)
257  {
258  perturbEnd = true;
259  }
260 
261 
262  if (perturbStart || perturbEnd)
263  {
264  const edge& e = surf1.edges()[edgeI];
265 
266  label v0 = surf1.meshPoints()[e[0]];
267  label v1 = surf1.meshPoints()[e[1]];
268 
269  vector eVec(points1[v1] - points1[v0]);
270  vector n = eVec/mag(eVec);
271 
272  if (perturbStart)
273  {
274  // Perturb with something (hopefully) larger than tolerance.
275  scalar t = 4.0*(rndGen.scalar01() - 0.5);
276  points1[v0] += t*surf1PointTol[e[0]]*n;
277 
278  const labelList& pEdges = surf1.pointEdges()[e[0]];
279 
280  forAll(pEdges, i)
281  {
282  affectedEdges[pEdges[i]] = true;
283  }
284  }
285  if (perturbEnd)
286  {
287  // Perturb with something larger than tolerance.
288  scalar t = 4.0*(rndGen.scalar01() - 0.5);
289  points1[v1] += t*surf1PointTol[e[1]]*n;
290 
291  const labelList& pEdges = surf1.pointEdges()[e[1]];
292 
293  forAll(pEdges, i)
294  {
295  affectedEdges[pEdges[i]] = true;
296  }
297  }
298  hasPerturbed = true;
299  }
300  }
301 
302  return hasPerturbed;
303 }
304 
305 
306 // Perturb single edge endpoint when perpendicular to face
307 bool Foam::edgeIntersections::rotatePerturb
308 (
309  const triSurface& surf1,
310  const scalarField& surf1PointTol, // surf1 tolerance per point
311  const label edgeI,
312 
313  Random& rndGen,
314  pointField& points1,
315  boolList& affectedEdges
316 ) const
317 {
318  const labelList& meshPoints = surf1.meshPoints();
319 
320  const labelList& edgeEnds = classification_[edgeI];
321 
322  bool hasPerturbed = false;
323 
324  forAll(edgeEnds, i)
325  {
326  if (edgeEnds[i] == 2)
327  {
328  const edge& e = surf1.edges()[edgeI];
329 
330  // Endpoint to modify. Choose either start or end.
331  label pointi = e[rndGen.sampleAB<label>(0, 2)];
332  // label pointi = e[0];
333 
334  // Generate random vector slightly larger than tolerance.
335  vector rndVec = rndGen.sample01<vector>() - vector(0.5, 0.5, 0.5);
336 
337  // Make sure rndVec only perp to edge
338  vector n(points1[meshPoints[e[1]]] - points1[meshPoints[e[0]]]);
339  scalar magN = mag(n) + vSmall;
340  n /= magN;
341 
342  rndVec -= n*(n & rndVec);
343 
344  // Normalize
345  rndVec /= mag(rndVec) + vSmall;
346 
347  // Scale to be moved by tolerance.
348  rndVec *= 0.01*magN;
349 
350  Pout<< "rotating: shifting endpoint " << meshPoints[pointi]
351  << " of edge:" << edgeI << " verts:"
352  << points1[meshPoints[e[0]]] << ' '
353  << points1[meshPoints[e[1]]]
354  << " by " << rndVec
355  << " tol:" << surf1PointTol[pointi] << endl;
356 
357  points1[meshPoints[pointi]] += rndVec;
358 
359  // Mark edges affected by change to point
360  const labelList& pEdges = surf1.pointEdges()[pointi];
361 
362  forAll(pEdges, i)
363  {
364  affectedEdges[pEdges[i]] = true;
365  }
366 
367  hasPerturbed = true;
368 
369  // Enough done for current edge; no need to test other intersections
370  // of this edge.
371  break;
372  }
373  }
374 
375  return hasPerturbed;
376 }
377 
378 
379 // Perturb edge when close to face
380 bool Foam::edgeIntersections::offsetPerturb
381 (
382  const triSurface& surf1,
383  const triSurface& surf2,
384  const label edgeI,
385 
386  Random& rndGen,
387  pointField& points1,
388  boolList& affectedEdges
389 ) const
390 {
391  const labelList& meshPoints = surf1.meshPoints();
392 
393  const edge& e = surf1.edges()[edgeI];
394 
395  const List<pointIndexHit>& hits = operator[](edgeI);
396 
397 
398  bool hasPerturbed = false;
399 
400  // For all hits on edge
401  forAll(hits, i)
402  {
403  const pointIndexHit& pHit = hits[i];
404 
405  // Classify point on face of surface2
406  label surf2Facei = pHit.index();
407 
408  const triSurface::FaceType& f2 = surf2.localFaces()[surf2Facei];
409  const pointField& surf2Pts = surf2.localPoints();
410 
411  const point ctr = f2.centre(surf2Pts);
412 
413  label nearType, nearLabel;
414 
415  f2.nearestPointClassify(pHit.hitPoint(), surf2Pts, nearType, nearLabel);
416 
417  if (nearType == triPointRef::POINT || nearType == triPointRef::EDGE)
418  {
419  // Shift edge towards tri centre
420  vector offset = 0.01*rndGen.scalar01()*(ctr - pHit.hitPoint());
421 
422  // shift e[0]
423  points1[meshPoints[e[0]]] += offset;
424 
425  // Mark edges affected by change to e0
426  const labelList& pEdges0 = surf1.pointEdges()[e[0]];
427 
428  forAll(pEdges0, i)
429  {
430  affectedEdges[pEdges0[i]] = true;
431  }
432 
433  // shift e[1]
434  points1[meshPoints[e[1]]] += offset;
435 
436  // Mark edges affected by change to e1
437  const labelList& pEdges1 = surf1.pointEdges()[e[1]];
438 
439  forAll(pEdges1, i)
440  {
441  affectedEdges[pEdges1[i]] = true;
442  }
443 
444  hasPerturbed = true;
445 
446  // No need to test any other hits on this edge
447  break;
448  }
449  }
450 
451  return hasPerturbed;
452 }
453 
454 
455 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
456 
457 // Construct null
459 :
460  List<List<pointIndexHit>>(),
461  classification_()
462 {}
463 
464 
465 // Construct from surface and tolerance
467 (
468  const triSurface& surf1,
469  const triSurfaceSearch& query2,
470  const scalarField& surf1PointTol
471 )
472 :
474  classification_(surf1.nEdges())
475 {
476  checkEdges(surf1);
477  checkEdges(query2.surface());
478 
479  // Current set of edges to test
480  labelList edgesToTest(surf1.nEdges());
481 
482  // Start off with all edges
483  forAll(edgesToTest, i)
484  {
485  edgesToTest[i] = i;
486  }
487 
488  // Determine intersections for edgesToTest
489  intersectEdges
490  (
491  surf1,
492  surf1.points(), // surf1 meshPoints (not localPoints!)
493  query2,
494  surf1PointTol,
495  edgesToTest
496  );
497 }
498 
499 
500 // Construct from components
502 (
503  const List<List<pointIndexHit>>& intersections,
504  const labelListList& classification
505 )
506 :
507  List<List<pointIndexHit>>(intersections),
508  classification_(classification)
509 {}
510 
511 
512 // * * * * * * * * * * * * * * * Static Functions * * * * * * * * * * * * * //
513 
515 {
516  const pointField& localPoints = surf.localPoints();
517  const labelListList& pointEdges = surf.pointEdges();
518  const edgeList& edges = surf.edges();
519 
520  scalarField minLen(localPoints.size());
521 
522  forAll(minLen, pointi)
523  {
524  const labelList& pEdges = pointEdges[pointi];
525 
526  scalar minDist = great;
527 
528  forAll(pEdges, i)
529  {
530  minDist = min(minDist, edges[pEdges[i]].mag(localPoints));
531  }
532 
533  minLen[pointi] = minDist;
534  }
535  return minLen;
536 }
537 
538 
539 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
540 
542 (
543  const label nIters,
544  const triSurface& surf1,
545  const triSurfaceSearch& query2,
546  const scalarField& surf1PointTol,
547  pointField& points1
548 )
549 {
550  const triSurface& surf2 = query2.surface();
551 
552  Random rndGen(356574);
553 
554  // Current set of edges to (re)test
555  labelList edgesToTest(surf1.nEdges());
556 
557  // Start off with all edges
558  forAll(edgesToTest, i)
559  {
560  edgesToTest[i] = i;
561  }
562 
563 
564  label iter = 0;
565 
566  for (; iter < nIters; iter++)
567  {
568  // Go through all edges to (re)test and perturb points if they are
569  // degenerate hits. Mark off edges that need to be recalculated.
570 
571  boolList affectedEdges(surf1.nEdges(), false);
572  label nShifted = 0;
573  label nRotated = 0;
574  label nOffset = 0;
575 
576  forAll(edgesToTest, i)
577  {
578  label edgeI = edgesToTest[i];
579 
580  // If edge not already marked for retesting
581  if (!affectedEdges[edgeI])
582  {
583  // 1. Check edges close to endpoint and perturb if necessary.
584 
585  bool shiftedEdgeEndPoints =
586  inlinePerturb
587  (
588  surf1,
589  surf1PointTol,
590  edgeI,
591  rndGen,
592  points1,
593  affectedEdges
594  );
595 
596  nShifted += (shiftedEdgeEndPoints ? 1 : 0);
597 
598  if (!shiftedEdgeEndPoints)
599  {
600  bool rotatedEdge =
601  rotatePerturb
602  (
603  surf1,
604  surf1PointTol,
605  edgeI,
606  rndGen,
607  points1,
608  affectedEdges
609  );
610 
611  nRotated += (rotatedEdge ? 1 : 0);
612 
613  if (!rotatedEdge)
614  {
615  // 2. we're sure now that the edge actually pierces the
616  // face. Now check the face for intersections close its
617  // points/edges
618 
619  bool offsetEdgePoints =
620  offsetPerturb
621  (
622  surf1,
623  surf2,
624  edgeI,
625  rndGen,
626  points1,
627  affectedEdges
628  );
629 
630  nOffset += (offsetEdgePoints ? 1 : 0);
631  }
632  }
633  }
634  }
635 
636  if (debug)
637  {
638  Pout<< "Edges to test : " << nl
639  << " total:" << edgesToTest.size() << nl
640  << " resolved by:" << nl
641  << " shifting : " << nShifted << nl
642  << " rotating : " << nRotated << nl
643  << " offsetting : " << nOffset << nl
644  << endl;
645  }
646 
647 
648  if (nShifted == 0 && nRotated == 0 && nOffset == 0)
649  {
650  // Nothing changed in current iteration. Current hit pattern is
651  // without any degenerates.
652  break;
653  }
654 
655  // Repack affected edges
656  labelList newEdgesToTest(surf1.nEdges());
657  label newEdgeI = 0;
658 
659  forAll(affectedEdges, edgeI)
660  {
661  if (affectedEdges[edgeI])
662  {
663  newEdgesToTest[newEdgeI++] = edgeI;
664  }
665  }
666  newEdgesToTest.setSize(newEdgeI);
667 
668  if (debug)
669  {
670  Pout<< "Edges to test:" << nl
671  << " was : " << edgesToTest.size() << nl
672  << " is : " << newEdgesToTest.size() << nl
673  << endl;
674  }
675 
676  // Transfer and test.
677  edgesToTest.transfer(newEdgesToTest);
678 
679  if (edgesToTest.empty())
680  {
681  FatalErrorInFunction << "oops" << abort(FatalError);
682  }
683 
684  // Re intersect moved edges.
685  intersectEdges
686  (
687  surf1,
688  points1, // surf1 meshPoints (not localPoints!)
689  query2,
690  surf1PointTol,
691  edgesToTest
692  );
693  }
694 
695  return iter;
696 }
697 
698 
700 (
701  const edgeIntersections& subInfo,
702  const labelList& edgeMap,
703  const labelList& faceMap,
704  const bool merge
705 )
706 {
707  forAll(subInfo, subI)
708  {
709  const List<pointIndexHit>& subHits = subInfo[subI];
710  const labelList& subClass = subInfo.classification()[subI];
711 
712  label edgeI = edgeMap[subI];
713  List<pointIndexHit>& intersections = operator[](edgeI);
714  labelList& intersectionTypes = classification_[edgeI];
715 
716  // Count unique hits. Assume edge can hit face only once
717  label sz = 0;
718  if (merge)
719  {
720  sz = intersections.size();
721  }
722 
723  label nNew = 0;
724  forAll(subHits, i)
725  {
726  const pointIndexHit& subHit = subHits[i];
727 
728  bool foundFace = false;
729  for (label interI = 0; interI < sz; interI++)
730  {
731  if (intersections[interI].index() == faceMap[subHit.index()])
732  {
733  foundFace = true;
734  break;
735  }
736  }
737  if (!foundFace)
738  {
739  nNew++;
740  }
741  }
742 
743 
744  intersections.setSize(sz+nNew);
745  intersectionTypes.setSize(sz+nNew);
746  nNew = sz;
747 
748  forAll(subHits, i)
749  {
750  const pointIndexHit& subHit = subHits[i];
751 
752  bool foundFace = false;
753  for (label interI = 0; interI < sz; interI++)
754  {
755  if (intersections[interI].index() == faceMap[subHit.index()])
756  {
757  foundFace = true;
758  break;
759  }
760  }
761 
762  if (!foundFace)
763  {
764  intersections[nNew] = pointIndexHit
765  (
766  subHit.hit(),
767  subHit.rawPoint(),
768  faceMap[subHit.index()]
769  );
770  intersectionTypes[nNew] = subClass[i];
771  nNew++;
772  }
773  }
774  }
775 }
776 
777 
778 // ************************************************************************* //
const labelListList & pointEdges() const
Return point-edge addressing.
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
error FatalError
T & operator[](const label)
Return element of UList.
Definition: UListI.H:167
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
Unit conversion functions.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
void merge(const edgeIntersections &, const labelList &edgeMap, const labelList &faceMap, const bool merge=true)
Merge (or override) edge intersection for a subset.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
scalar minDist(const List< pointIndexHit > &hitList)
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
Random rndGen(label(0))
const Vector< Cmpt > & centre(const Foam::List< Vector< Cmpt >> &) const
Return *this (used for point which is a typedef to Vector<scalar>.
Definition: VectorI.H:116
Helper class to search on triSurface.
static scalarField minEdgeLength(const triSurface &surf)
Calculate min edge length for every surface point.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dimensionedScalar cos(const dimensionedScalar &ds)
List< edge > edgeList
Definition: edgeList.H:38
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
bool hit() const
Is there a hit.
const Field< PointType > & points() const
Return reference to global points.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
List< label > labelList
A List of labels.
Definition: labelList.H:56
const labelListList & classification() const
For every intersection the classification status.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Random number generator.
Definition: Random.H:57
const Point & rawPoint() const
Return point with no checking.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
label nEdges() const
Return number of edges in patch.
static const char nl
Definition: Ostream.H:265
defineTypeNameAndDebug(combustionModel, 0)
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.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const triSurface & surface() const
Return reference to the surface.
void setSize(const label)
Reset size of List.
Definition: List.C:281
vector point
Point is a vector.
Definition: point.H:41
#define WarningInFunction
Report a warning using Foam::Warning.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
Triangulated surface description with patch information.
Definition: triSurface.H:66
edgeIntersections()
Construct null.
label index() const
Return index.
Namespace for OpenFOAM.