triSurfaceTools.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 | Copyright (C) 2011-2015 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 "triSurfaceTools.H"
27 
28 #include "triSurface.H"
29 #include "OFstream.H"
30 #include "mergePoints.H"
31 #include "polyMesh.H"
32 #include "plane.H"
33 #include "geompack.H"
34 
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 /*
45  Refine by splitting all three edges of triangle ('red' refinement).
46  Neighbouring triangles (which are not marked for refinement get split
47  in half ('green') refinement. (R. Verfuerth, "A review of a posteriori
48  error estimation and adaptive mesh refinement techniques",
49  Wiley-Teubner, 1996)
50 */
51 
52 // FaceI gets refined ('red'). Propagate edge refinement.
54 (
55  const triSurface& surf,
56  const label faceI,
57  List<refineType>& refine
58 )
59 {
60  if (refine[faceI] == RED)
61  {
62  // Already marked for refinement. Do nothing.
63  }
64  else
65  {
66  // Not marked or marked for 'green' refinement. Refine.
67  refine[faceI] = RED;
68 
69  const labelList& myNeighbours = surf.faceFaces()[faceI];
70 
71  forAll(myNeighbours, myNeighbourI)
72  {
73  label neighbourFaceI = myNeighbours[myNeighbourI];
74 
75  if (refine[neighbourFaceI] == GREEN)
76  {
77  // Change to red refinement and propagate
78  calcRefineStatus(surf, neighbourFaceI, refine);
79  }
80  else if (refine[neighbourFaceI] == NONE)
81  {
82  refine[neighbourFaceI] = GREEN;
83  }
84  }
85  }
86 }
87 
88 
89 // Split faceI along edgeI at position newPointI
91 (
92  const triSurface& surf,
93  const label faceI,
94  const label edgeI,
95  const label newPointI,
96  DynamicList<labelledTri>& newFaces
97 )
98 {
99  const labelledTri& f = surf.localFaces()[faceI];
100  const edge& e = surf.edges()[edgeI];
101 
102  // Find index of edge in face.
103 
104  label fp0 = findIndex(f, e[0]);
105  label fp1 = f.fcIndex(fp0);
106  label fp2 = f.fcIndex(fp1);
107 
108  if (f[fp1] == e[1])
109  {
110  // Edge oriented like face
111  newFaces.append
112  (
114  (
115  f[fp0],
116  newPointI,
117  f[fp2],
118  f.region()
119  )
120  );
121  newFaces.append
122  (
124  (
125  newPointI,
126  f[fp1],
127  f[fp2],
128  f.region()
129  )
130  );
131  }
132  else
133  {
134  newFaces.append
135  (
137  (
138  f[fp2],
139  newPointI,
140  f[fp1],
141  f.region()
142  )
143  );
144  newFaces.append
145  (
147  (
148  newPointI,
149  f[fp0],
150  f[fp1],
151  f.region()
152  )
153  );
154  }
155 }
156 
157 
158 // Refine all triangles marked for refinement.
160 (
161  const triSurface& surf,
162  const List<refineType>& refineStatus
163 )
164 {
165  // Storage for new points. (start after old points)
166  DynamicList<point> newPoints(surf.nPoints());
167  forAll(surf.localPoints(), pointI)
168  {
169  newPoints.append(surf.localPoints()[pointI]);
170  }
171  label newVertI = surf.nPoints();
172 
173  // Storage for new faces
174  DynamicList<labelledTri> newFaces(surf.size());
175 
176 
177  // Point index for midpoint on edge
178  labelList edgeMid(surf.nEdges(), -1);
179 
180  forAll(refineStatus, faceI)
181  {
182  if (refineStatus[faceI] == RED)
183  {
184  // Create new vertices on all edges to be refined.
185  const labelList& fEdges = surf.faceEdges()[faceI];
186 
187  forAll(fEdges, i)
188  {
189  label edgeI = fEdges[i];
190 
191  if (edgeMid[edgeI] == -1)
192  {
193  const edge& e = surf.edges()[edgeI];
194 
195  // Create new point on mid of edge
196  newPoints.append
197  (
198  0.5
199  * (
200  surf.localPoints()[e.start()]
201  + surf.localPoints()[e.end()]
202  )
203  );
204  edgeMid[edgeI] = newVertI++;
205  }
206  }
207 
208  // Now we have new mid edge vertices for all edges on face
209  // so create triangles for RED rerfinement.
210 
211  const edgeList& edges = surf.edges();
212 
213  // Corner triangles
214  newFaces.append
215  (
217  (
218  edgeMid[fEdges[0]],
219  edges[fEdges[0]].commonVertex(edges[fEdges[1]]),
220  edgeMid[fEdges[1]],
221  surf[faceI].region()
222  )
223  );
224 
225  newFaces.append
226  (
228  (
229  edgeMid[fEdges[1]],
230  edges[fEdges[1]].commonVertex(edges[fEdges[2]]),
231  edgeMid[fEdges[2]],
232  surf[faceI].region()
233  )
234  );
235 
236  newFaces.append
237  (
239  (
240  edgeMid[fEdges[2]],
241  edges[fEdges[2]].commonVertex(edges[fEdges[0]]),
242  edgeMid[fEdges[0]],
243  surf[faceI].region()
244  )
245  );
246 
247  // Inner triangle
248  newFaces.append
249  (
251  (
252  edgeMid[fEdges[0]],
253  edgeMid[fEdges[1]],
254  edgeMid[fEdges[2]],
255  surf[faceI].region()
256  )
257  );
258 
259 
260  // Create triangles for GREEN refinement.
261  forAll(fEdges, i)
262  {
263  const label edgeI = fEdges[i];
264 
265  label otherFaceI = otherFace(surf, faceI, edgeI);
266 
267  if ((otherFaceI != -1) && (refineStatus[otherFaceI] == GREEN))
268  {
269  greenRefine
270  (
271  surf,
272  otherFaceI,
273  edgeI,
274  edgeMid[edgeI],
275  newFaces
276  );
277  }
278  }
279  }
280  }
281 
282  // Copy unmarked triangles since keep original vertex numbering.
283  forAll(refineStatus, faceI)
284  {
285  if (refineStatus[faceI] == NONE)
286  {
287  newFaces.append(surf.localFaces()[faceI]);
288  }
289  }
290 
291  newFaces.shrink();
292  newPoints.shrink();
293 
294 
295  // Transfer DynamicLists to straight ones.
297  allPoints.transfer(newPoints);
298  newPoints.clear();
299 
300  return triSurface(newFaces, surf.patches(), allPoints, true);
301 }
302 
303 
304 // Angle between two neighbouring triangles,
305 // angle between shared-edge end points and left and right face end points
307 (
308  const point& pStart,
309  const point& pEnd,
310  const point& pLeft,
311  const point& pRight
312 )
313 {
314  const vector common(pEnd - pStart);
315  const vector base0(pLeft - pStart);
316  const vector base1(pRight - pStart);
317 
318  vector n0(common ^ base0);
319  n0 /= Foam::mag(n0);
320 
321  vector n1(base1 ^ common);
322  n1 /= Foam::mag(n1);
323 
324  return n0 & n1;
325 }
326 
327 
328 // Protect edges around vertex from collapsing.
329 // Moving vertI to a different position
330 // - affects obviously the cost of the faces using it
331 // - but also their neighbours since the edge between the faces changes
333 (
334  const triSurface& surf,
335  const label vertI,
336  labelList& faceStatus
337 )
338 {
339 // const labelList& myFaces = surf.pointFaces()[vertI];
340 // forAll(myFaces, i)
341 // {
342 // label faceI = myFaces[i];
343 //
344 // if ((faceStatus[faceI] == ANYEDGE) || (faceStatus[faceI] >= 0))
345 // {
346 // faceStatus[faceI] = NOEDGE;
347 // }
348 // }
349 
350  const labelList& myEdges = surf.pointEdges()[vertI];
351  forAll(myEdges, i)
352  {
353  const labelList& myFaces = surf.edgeFaces()[myEdges[i]];
354 
355  forAll(myFaces, myFaceI)
356  {
357  label faceI = myFaces[myFaceI];
358 
359  if ((faceStatus[faceI] == ANYEDGE) || (faceStatus[faceI] >= 0))
360  {
361  faceStatus[faceI] = NOEDGE;
362  }
363  }
364  }
365 }
366 
367 
368 //
369 // Edge collapse helper functions
370 //
371 
372 // Get all faces that will get collapsed if edgeI collapses.
374 (
375  const triSurface& surf,
376  label edgeI
377 )
378 {
379  const edge& e = surf.edges()[edgeI];
380  label v1 = e.start();
381  label v2 = e.end();
382 
383  // Faces using edge will certainly get collapsed.
384  const labelList& myFaces = surf.edgeFaces()[edgeI];
385 
386  labelHashSet facesToBeCollapsed(2*myFaces.size());
387 
388  forAll(myFaces, myFaceI)
389  {
390  facesToBeCollapsed.insert(myFaces[myFaceI]);
391  }
392 
393  // From faces using v1 check if they share an edge with faces
394  // using v2.
395  // - share edge: are part of 'splay' tree and will collapse if edge
396  // collapses
397  const labelList& v1Faces = surf.pointFaces()[v1];
398 
399  forAll(v1Faces, v1FaceI)
400  {
401  label face1I = v1Faces[v1FaceI];
402 
403  label otherEdgeI = oppositeEdge(surf, face1I, v1);
404 
405  // Step across edge to other face
406  label face2I = otherFace(surf, face1I, otherEdgeI);
407 
408  if (face2I != -1)
409  {
410  // found face on other side of edge. Now check if uses v2.
411  if (oppositeVertex(surf, face2I, otherEdgeI) == v2)
412  {
413  // triangles face1I and face2I form splay tree and will
414  // collapse.
415  facesToBeCollapsed.insert(face1I);
416  facesToBeCollapsed.insert(face2I);
417  }
418  }
419  }
420 
421  return facesToBeCollapsed;
422 }
423 
424 
425 // Return value of faceUsed for faces using vertI
427 (
428  const triSurface& surf,
429  const labelHashSet& faceUsed,
430  const label vertI
431 )
432 {
433  const labelList& myFaces = surf.pointFaces()[vertI];
434 
435  forAll(myFaces, myFaceI)
436  {
437  label face1I = myFaces[myFaceI];
438 
439  if (faceUsed.found(face1I))
440  {
441  return face1I;
442  }
443  }
444  return -1;
445 }
446 
447 
448 // Calculate new edge-face addressing resulting from edge collapse
450 (
451  const triSurface& surf,
452  const label edgeI,
453  const labelHashSet& collapsedFaces,
454  HashTable<label, label, Hash<label> >& edgeToEdge,
455  HashTable<label, label, Hash<label> >& edgeToFace
456 )
457 {
458  const edge& e = surf.edges()[edgeI];
459  label v1 = e.start();
460  label v2 = e.end();
461 
462  const labelList& v1Faces = surf.pointFaces()[v1];
463  const labelList& v2Faces = surf.pointFaces()[v2];
464 
465  // Mark all (non collapsed) faces using v2
466  labelHashSet v2FacesHash(v2Faces.size());
467 
468  forAll(v2Faces, v2FaceI)
469  {
470  if (!collapsedFaces.found(v2Faces[v2FaceI]))
471  {
472  v2FacesHash.insert(v2Faces[v2FaceI]);
473  }
474  }
475 
476 
477  forAll(v1Faces, v1FaceI)
478  {
479  label face1I = v1Faces[v1FaceI];
480 
481  if (collapsedFaces.found(face1I))
482  {
483  continue;
484  }
485 
486  //
487  // Check if face1I uses a vertex connected to a face using v2
488  //
489 
490  label vert1I = -1;
491  label vert2I = -1;
492  otherVertices
493  (
494  surf,
495  face1I,
496  v1,
497  vert1I,
498  vert2I
499  );
500  //Pout<< "Face:" << surf.localFaces()[face1I] << " other vertices:"
501  // << vert1I << ' ' << vert2I << endl;
502 
503  // Check vert1, vert2 for usage by v2Face.
504 
505  label commonVert = vert1I;
506  label face2I = vertexUsesFace(surf, v2FacesHash, commonVert);
507  if (face2I == -1)
508  {
509  commonVert = vert2I;
510  face2I = vertexUsesFace(surf, v2FacesHash, commonVert);
511  }
512 
513  if (face2I != -1)
514  {
515  // Found one: commonVert is used by both face1I and face2I
516  label edge1I = getEdge(surf, v1, commonVert);
517  label edge2I = getEdge(surf, v2, commonVert);
518 
519  edgeToEdge.insert(edge1I, edge2I);
520  edgeToEdge.insert(edge2I, edge1I);
521 
522  edgeToFace.insert(edge1I, face2I);
523  edgeToFace.insert(edge2I, face1I);
524  }
525  }
526 }
527 
528 
529 // Calculates (cos of) angle across edgeI of faceI,
530 // taking into account updated addressing (resulting from edge collapse)
532 (
533  const triSurface& surf,
534  const label v1,
535  const point& pt,
536  const labelHashSet& collapsedFaces,
537  const HashTable<label, label, Hash<label> >& edgeToEdge,
538  const HashTable<label, label, Hash<label> >& edgeToFace,
539  const label faceI,
540  const label edgeI
541 )
542 {
543  const pointField& localPoints = surf.localPoints();
544 
545  label A = surf.edges()[edgeI].start();
546  label B = surf.edges()[edgeI].end();
547  label C = oppositeVertex(surf, faceI, edgeI);
548 
549  label D = -1;
550 
551  label face2I = -1;
552 
553  if (edgeToEdge.found(edgeI))
554  {
555  // Use merged addressing
556  label edge2I = edgeToEdge[edgeI];
557  face2I = edgeToFace[edgeI];
558 
559  D = oppositeVertex(surf, face2I, edge2I);
560  }
561  else
562  {
563  // Use normal edge-face addressing
564  face2I = otherFace(surf, faceI, edgeI);
565 
566  if ((face2I != -1) && !collapsedFaces.found(face2I))
567  {
568  D = oppositeVertex(surf, face2I, edgeI);
569  }
570  }
571 
572  scalar cosAngle = 1;
573 
574  if (D != -1)
575  {
576  if (A == v1)
577  {
578  cosAngle = faceCosAngle
579  (
580  pt,
581  localPoints[B],
582  localPoints[C],
583  localPoints[D]
584  );
585  }
586  else if (B == v1)
587  {
588  cosAngle = faceCosAngle
589  (
590  localPoints[A],
591  pt,
592  localPoints[C],
593  localPoints[D]
594  );
595  }
596  else if (C == v1)
597  {
598  cosAngle = faceCosAngle
599  (
600  localPoints[A],
601  localPoints[B],
602  pt,
603  localPoints[D]
604  );
605  }
606  else if (D == v1)
607  {
608  cosAngle = faceCosAngle
609  (
610  localPoints[A],
611  localPoints[B],
612  localPoints[C],
613  pt
614  );
615  }
616  else
617  {
619  << "face " << faceI << " does not use vertex "
620  << v1 << " of collapsed edge" << abort(FatalError);
621  }
622  }
623  return cosAngle;
624 }
625 
626 
627 //- Calculate minimum (cos of) edge angle using addressing from collapsing
628 // edge to v1 at pt.
630 (
631  const triSurface& surf,
632  const label v1,
633  const point& pt,
634  const labelHashSet& collapsedFaces,
635  const HashTable<label, label, Hash<label> >& edgeToEdge,
636  const HashTable<label, label, Hash<label> >& edgeToFace
637 )
638 {
639  const labelList& v1Faces = surf.pointFaces()[v1];
640 
641  scalar minCos = 1;
642 
643  forAll(v1Faces, v1FaceI)
644  {
645  label faceI = v1Faces[v1FaceI];
646 
647  if (collapsedFaces.found(faceI))
648  {
649  continue;
650  }
651 
652  const labelList& myEdges = surf.faceEdges()[faceI];
653 
654  forAll(myEdges, myEdgeI)
655  {
656  label edgeI = myEdges[myEdgeI];
657 
658  minCos =
659  min
660  (
661  minCos,
662  edgeCosAngle
663  (
664  surf,
665  v1,
666  pt,
667  collapsedFaces,
668  edgeToEdge,
669  edgeToFace,
670  faceI,
671  edgeI
672  )
673  );
674  }
675  }
676 
677  return minCos;
678 }
679 
680 
681 // Calculate max dihedral angle after collapsing edge to v1 (at pt).
682 // Note that all edges of all faces using v1 are affected.
684 (
685  const triSurface& surf,
686  const label v1,
687  const point& pt,
688  const labelHashSet& collapsedFaces,
689  const HashTable<label, label, Hash<label> >& edgeToEdge,
690  const HashTable<label, label, Hash<label> >& edgeToFace,
691  const scalar minCos
692 )
693 {
694  const labelList& v1Faces = surf.pointFaces()[v1];
695 
696  forAll(v1Faces, v1FaceI)
697  {
698  label faceI = v1Faces[v1FaceI];
699 
700  if (collapsedFaces.found(faceI))
701  {
702  continue;
703  }
704 
705  const labelList& myEdges = surf.faceEdges()[faceI];
706 
707  forAll(myEdges, myEdgeI)
708  {
709  label edgeI = myEdges[myEdgeI];
710 
711  if
712  (
713  edgeCosAngle
714  (
715  surf,
716  v1,
717  pt,
718  collapsedFaces,
719  edgeToEdge,
720  edgeToFace,
721  faceI,
722  edgeI
723  )
724  < minCos
725  )
726  {
727  return true;
728  }
729  }
730  }
731 
732  return false;
733 }
734 
735 
738 // a vertex.
739 //bool Foam::triSurfaceTools::collapseCreatesDuplicates
740 //(
741 // const triSurface& surf,
742 // const label edgeI,
743 // const labelHashSet& collapsedFaces
744 //)
745 //{
746 //Pout<< "duplicate : edgeI:" << edgeI << surf.edges()[edgeI]
747 // << " collapsedFaces:" << collapsedFaces.toc() << endl;
748 //
749 // // Mark neighbours of faces to be collapsed, i.e. get the first layer
750 // // of triangles outside collapsedFaces.
751 // // neighbours actually contains the
752 // // edge with which triangle connects to collapsedFaces.
753 //
754 // HashTable<label, label, Hash<label> > neighbours;
755 //
756 // labelList collapsed = collapsedFaces.toc();
757 //
758 // forAll(collapsed, collapseI)
759 // {
760 // const label faceI = collapsed[collapseI];
761 //
762 // const labelList& myEdges = surf.faceEdges()[faceI];
763 //
764 // Pout<< "collapsing faceI:" << faceI << " uses edges:" << myEdges
765 // << endl;
766 //
767 // forAll(myEdges, myEdgeI)
768 // {
769 // const labelList& myFaces = surf.edgeFaces()[myEdges[myEdgeI]];
770 //
771 // Pout<< "Edge " << myEdges[myEdgeI] << " is used by faces "
772 // << myFaces << endl;
773 //
774 // if ((myEdges[myEdgeI] != edgeI) && (myFaces.size() == 2))
775 // {
776 // // Get the other face
777 // label neighbourFaceI = myFaces[0];
778 //
779 // if (neighbourFaceI == faceI)
780 // {
781 // neighbourFaceI = myFaces[1];
782 // }
783 //
784 // // Is 'outside' face if not in collapsedFaces itself
785 // if (!collapsedFaces.found(neighbourFaceI))
786 // {
787 // neighbours.insert(neighbourFaceI, myEdges[myEdgeI]);
788 // }
789 // }
790 // }
791 // }
792 //
793 // // Now neighbours contains first layer of triangles outside of
794 // // collapseFaces
795 // // There should be
796 // // -two if edgeI is a boundary edge
797 // // since the outside 'edge' of collapseFaces should
798 // // form a triangle and the face connected to edgeI is not inserted.
799 // // -four if edgeI is not a boundary edge since then the outside edge forms
800 // // a diamond.
801 //
802 // // Check if any of the faces in neighbours share an edge. (n^2)
803 //
804 // labelList neighbourList = neighbours.toc();
805 //
806 //Pout<< "edgeI:" << edgeI << " neighbourList:" << neighbourList << endl;
807 //
808 //
809 // forAll(neighbourList, i)
810 // {
811 // const labelList& faceIEdges = surf.faceEdges()[neighbourList[i]];
812 //
813 // for (label j = i+1; j < neighbourList.size(); i++)
814 // {
815 // const labelList& faceJEdges = surf.faceEdges()[neighbourList[j]];
816 //
817 // // Check if faceI and faceJ share an edge
818 // forAll(faceIEdges, fI)
819 // {
820 // forAll(faceJEdges, fJ)
821 // {
822 // Pout<< " comparing " << faceIEdges[fI] << " to "
823 // << faceJEdges[fJ] << endl;
824 //
825 // if (faceIEdges[fI] == faceJEdges[fJ])
826 // {
827 // return true;
828 // }
829 // }
830 // }
831 // }
832 // }
833 // Pout<< "Found no match. Returning false" << endl;
834 // return false;
835 //}
836 
837 
838 // Finds the triangle edge cut by the plane between a point inside / on edge
839 // of a triangle and a point outside. Returns:
840 // - cut through edge/point
841 // - miss()
843 (
844  const triSurface& s,
845  const label triI,
846  const label excludeEdgeI,
847  const label excludePointI,
848 
849  const point& triPoint,
850  const plane& cutPlane,
851  const point& toPoint
852 )
853 {
854  const pointField& points = s.points();
855  const labelledTri& f = s[triI];
856  const labelList& fEdges = s.faceEdges()[triI];
857 
858  // Get normal distance to planeN
860 
861  scalar norm = 0;
862  forAll(d, fp)
863  {
864  d[fp] = (points[f[fp]]-cutPlane.refPoint()) & cutPlane.normal();
865  norm += mag(d[fp]);
866  }
867 
868  // Normalise and truncate
869  forAll(d, i)
870  {
871  d[i] /= norm;
872 
873  if (mag(d[i]) < 1e-6)
874  {
875  d[i] = 0.0;
876  }
877  }
878 
879  // Return information
881 
882  if (excludePointI != -1)
883  {
884  // Excluded point. Test only opposite edge.
885 
886  label fp0 = findIndex(s.localFaces()[triI], excludePointI);
887 
888  if (fp0 == -1)
889  {
891  << " localF:" << s.localFaces()[triI] << abort(FatalError);
892  }
893 
894  label fp1 = f.fcIndex(fp0);
895  label fp2 = f.fcIndex(fp1);
896 
897 
898  if (d[fp1] == 0.0)
899  {
900  cut.setHit();
901  cut.setPoint(points[f[fp1]]);
902  cut.elementType() = triPointRef::POINT;
903  cut.setIndex(s.localFaces()[triI][fp1]);
904  }
905  else if (d[fp2] == 0.0)
906  {
907  cut.setHit();
908  cut.setPoint(points[f[fp2]]);
909  cut.elementType() = triPointRef::POINT;
910  cut.setIndex(s.localFaces()[triI][fp2]);
911  }
912  else if
913  (
914  (d[fp1] < 0 && d[fp2] < 0)
915  || (d[fp1] > 0 && d[fp2] > 0)
916  )
917  {
918  // Both same sign. Not crossing edge at all.
919  // cut already set to miss().
920  }
921  else
922  {
923  cut.setHit();
924  cut.setPoint
925  (
926  (d[fp2]*points[f[fp1]] - d[fp1]*points[f[fp2]])
927  / (d[fp2] - d[fp1])
928  );
929  cut.elementType() = triPointRef::EDGE;
930  cut.setIndex(fEdges[fp1]);
931  }
932  }
933  else
934  {
935  // Find the two intersections
937  label interI = 0;
938 
939  forAll(f, fp0)
940  {
941  label fp1 = f.fcIndex(fp0);
942 
943  if (d[fp0] == 0)
944  {
945  if (interI >= 2)
946  {
948  << "problem : triangle has three intersections." << nl
949  << "triangle:" << f.tri(points)
950  << " d:" << d << abort(FatalError);
951  }
952  inters[interI].setHit();
953  inters[interI].setPoint(points[f[fp0]]);
954  inters[interI].elementType() = triPointRef::POINT;
955  inters[interI].setIndex(s.localFaces()[triI][fp0]);
956  interI++;
957  }
958  else if
959  (
960  (d[fp0] < 0 && d[fp1] > 0)
961  || (d[fp0] > 0 && d[fp1] < 0)
962  )
963  {
964  if (interI >= 2)
965  {
967  << "problem : triangle has three intersections." << nl
968  << "triangle:" << f.tri(points)
969  << " d:" << d << abort(FatalError);
970  }
971  inters[interI].setHit();
972  inters[interI].setPoint
973  (
974  (d[fp0]*points[f[fp1]] - d[fp1]*points[f[fp0]])
975  / (d[fp0] - d[fp1])
976  );
977  inters[interI].elementType() = triPointRef::EDGE;
978  inters[interI].setIndex(fEdges[fp0]);
979  interI++;
980  }
981  }
982 
983 
984  if (interI == 0)
985  {
986  // Return miss
987  }
988  else if (interI == 1)
989  {
990  // Only one intersection. Should not happen!
991  cut = inters[0];
992  }
993  else if (interI == 2)
994  {
995  // Handle excludeEdgeI
996  if
997  (
998  inters[0].elementType() == triPointRef::EDGE
999  && inters[0].index() == excludeEdgeI
1000  )
1001  {
1002  cut = inters[1];
1003  }
1004  else if
1005  (
1006  inters[1].elementType() == triPointRef::EDGE
1007  && inters[1].index() == excludeEdgeI
1008  )
1009  {
1010  cut = inters[0];
1011  }
1012  else
1013  {
1014  // Two cuts. Find nearest.
1015  if
1016  (
1017  magSqr(inters[0].rawPoint() - toPoint)
1018  < magSqr(inters[1].rawPoint() - toPoint)
1019  )
1020  {
1021  cut = inters[0];
1022  }
1023  else
1024  {
1025  cut = inters[1];
1026  }
1027  }
1028  }
1029  }
1030  return cut;
1031 }
1032 
1033 
1034 // 'Snap' point on to endPoint.
1037  const triSurface& s,
1038  const surfaceLocation& end,
1039  surfaceLocation& current
1040 )
1041 {
1042  if (end.elementType() == triPointRef::NONE)
1043  {
1044  if (current.elementType() == triPointRef::NONE)
1045  {
1046  // endpoint on triangle; current on triangle
1047  if (current.index() == end.index())
1048  {
1049  //if (debug)
1050  //{
1051  // Pout<< "snapToEnd : snapping:" << current << " onto:"
1052  // << end << endl;
1053  //}
1054  current = end;
1055  current.setHit();
1056  }
1057  }
1058  // No need to handle current on edge/point since tracking handles this.
1059  }
1060  else if (end.elementType() == triPointRef::EDGE)
1061  {
1062  if (current.elementType() == triPointRef::NONE)
1063  {
1064  // endpoint on edge; current on triangle
1065  const labelList& fEdges = s.faceEdges()[current.index()];
1066 
1067  if (findIndex(fEdges, end.index()) != -1)
1068  {
1069  //if (debug)
1070  //{
1071  // Pout<< "snapToEnd : snapping:" << current << " onto:"
1072  // << end << endl;
1073  //}
1074  current = end;
1075  current.setHit();
1076  }
1077  }
1078  else if (current.elementType() == triPointRef::EDGE)
1079  {
1080  // endpoint on edge; current on edge
1081  if (current.index() == end.index())
1082  {
1083  //if (debug)
1084  //{
1085  // Pout<< "snapToEnd : snapping:" << current << " onto:"
1086  // << end << endl;
1087  //}
1088  current = end;
1089  current.setHit();
1090  }
1091  }
1092  else
1093  {
1094  // endpoint on edge; current on point
1095  const edge& e = s.edges()[end.index()];
1096 
1097  if (current.index() == e[0] || current.index() == e[1])
1098  {
1099  //if (debug)
1100  //{
1101  // Pout<< "snapToEnd : snapping:" << current << " onto:"
1102  // << end << endl;
1103  //}
1104  current = end;
1105  current.setHit();
1106  }
1107  }
1108  }
1109  else // end.elementType() == POINT
1110  {
1111  if (current.elementType() == triPointRef::NONE)
1112  {
1113  // endpoint on point; current on triangle
1114  const triSurface::FaceType& f = s.localFaces()[current.index()];
1115 
1116  if (findIndex(f, end.index()) != -1)
1117  {
1118  //if (debug)
1119  //{
1120  // Pout<< "snapToEnd : snapping:" << current << " onto:"
1121  // << end << endl;
1122  //}
1123  current = end;
1124  current.setHit();
1125  }
1126  }
1127  else if (current.elementType() == triPointRef::EDGE)
1128  {
1129  // endpoint on point; current on edge
1130  const edge& e = s.edges()[current.index()];
1131 
1132  if (end.index() == e[0] || end.index() == e[1])
1133  {
1134  //if (debug)
1135  //{
1136  // Pout<< "snapToEnd : snapping:" << current << " onto:"
1137  // << end << endl;
1138  //}
1139  current = end;
1140  current.setHit();
1141  }
1142  }
1143  else
1144  {
1145  // endpoint on point; current on point
1146  if (current.index() == end.index())
1147  {
1148  //if (debug)
1149  //{
1150  // Pout<< "snapToEnd : snapping:" << current << " onto:"
1151  // << end << endl;
1152  //}
1153  current = end;
1154  current.setHit();
1155  }
1156  }
1157  }
1158 }
1159 
1160 
1161 // Start:
1162 // - location
1163 // - element type (triangle/edge/point) and index
1164 // - triangle to exclude
1167  const triSurface& s,
1168  const labelList& eFaces,
1169  const surfaceLocation& start,
1170  const label excludeEdgeI,
1171  const label excludePointI,
1172  const surfaceLocation& end,
1173  const plane& cutPlane
1174 )
1175 {
1176  surfaceLocation nearest;
1177 
1178  scalar minDistSqr = Foam::sqr(GREAT);
1179 
1180  forAll(eFaces, i)
1181  {
1182  label triI = eFaces[i];
1183 
1184  // Make sure we don't revisit previous face
1185  if (triI != start.triangle())
1186  {
1187  if (end.elementType() == triPointRef::NONE && end.index() == triI)
1188  {
1189  // Endpoint is in this triangle. Jump there.
1190  nearest = end;
1191  nearest.setHit();
1192  nearest.triangle() = triI;
1193  break;
1194  }
1195  else
1196  {
1197  // Which edge is cut.
1198 
1199  surfaceLocation cutInfo = cutEdge
1200  (
1201  s,
1202  triI,
1203  excludeEdgeI, // excludeEdgeI
1204  excludePointI, // excludePointI
1205  start.rawPoint(),
1206  cutPlane,
1207  end.rawPoint()
1208  );
1209 
1210  // If crossing an edge we expect next edge to be cut.
1211  if (excludeEdgeI != -1 && !cutInfo.hit())
1212  {
1214  << "Triangle:" << triI
1215  << " excludeEdge:" << excludeEdgeI
1216  << " point:" << start.rawPoint()
1217  << " plane:" << cutPlane
1218  << " . No intersection!" << abort(FatalError);
1219  }
1220 
1221  if (cutInfo.hit())
1222  {
1223  scalar distSqr = magSqr(cutInfo.rawPoint()-end.rawPoint());
1224 
1225  if (distSqr < minDistSqr)
1226  {
1227  minDistSqr = distSqr;
1228  nearest = cutInfo;
1229  nearest.triangle() = triI;
1230  nearest.setMiss();
1231  }
1232  }
1233  }
1234  }
1235  }
1236 
1237  if (nearest.triangle() == -1)
1238  {
1239  // Did not move from edge. Give warning? Return something special?
1240  // For now responsability of caller to make sure that nothing has
1241  // moved.
1242  }
1243 
1244  return nearest;
1245 }
1246 
1247 
1248 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1249 
1250 
1251 // Write pointField to file
1254  const fileName& fName,
1255  const pointField& pts
1256 )
1257 {
1258  OFstream outFile(fName);
1259 
1260  forAll(pts, pointI)
1261  {
1262  const point& pt = pts[pointI];
1263 
1264  outFile<< "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
1265  }
1266  Pout<< "Written " << pts.size() << " vertices to file " << fName << endl;
1267 }
1268 
1269 
1270 // Write vertex subset to OBJ format file
1273  const triSurface& surf,
1274  const fileName& fName,
1275  const boolList& markedVerts
1276 )
1277 {
1278  OFstream outFile(fName);
1279 
1280  label nVerts = 0;
1281  forAll(markedVerts, vertI)
1282  {
1283  if (markedVerts[vertI])
1284  {
1285  const point& pt = surf.localPoints()[vertI];
1286 
1287  outFile<< "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
1288 
1289  nVerts++;
1290  }
1291  }
1292  Pout<< "Written " << nVerts << " vertices to file " << fName << endl;
1293 }
1294 
1295 
1296 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1297 // Addressing helper functions:
1298 
1299 
1300 // Get all triangles using vertices of edge
1303  const triSurface& surf,
1304  const label edgeI,
1305  labelList& edgeTris
1306 )
1307 {
1308  const edge& e = surf.edges()[edgeI];
1309  const labelList& myFaces = surf.edgeFaces()[edgeI];
1310 
1311  label face1I = myFaces[0];
1312  label face2I = -1;
1313  if (myFaces.size() == 2)
1314  {
1315  face2I = myFaces[1];
1316  }
1317 
1318  const labelList& startFaces = surf.pointFaces()[e.start()];
1319  const labelList& endFaces = surf.pointFaces()[e.end()];
1320 
1321  // Number of triangles is sum of pointfaces - common faces
1322  // (= faces using edge)
1323  edgeTris.setSize(startFaces.size() + endFaces.size() - myFaces.size());
1324 
1325  label nTris = 0;
1326  forAll(startFaces, startFaceI)
1327  {
1328  edgeTris[nTris++] = startFaces[startFaceI];
1329  }
1330 
1331  forAll(endFaces, endFaceI)
1332  {
1333  label faceI = endFaces[endFaceI];
1334 
1335  if ((faceI != face1I) && (faceI != face2I))
1336  {
1337  edgeTris[nTris++] = faceI;
1338  }
1339  }
1340 }
1341 
1342 
1343 // Get all vertices connected to vertices of edge
1346  const triSurface& surf,
1347  const edge& e
1348 )
1349 {
1350  const edgeList& edges = surf.edges();
1351 
1352  label v1 = e.start();
1353  label v2 = e.end();
1354 
1355  // Get all vertices connected to v1 or v2 through an edge
1356  labelHashSet vertexNeighbours;
1357 
1358  const labelList& v1Edges = surf.pointEdges()[v1];
1359 
1360  forAll(v1Edges, v1EdgeI)
1361  {
1362  const edge& e = edges[v1Edges[v1EdgeI]];
1363  vertexNeighbours.insert(e.otherVertex(v1));
1364  }
1365 
1366  const labelList& v2Edges = surf.pointEdges()[v2];
1367 
1368  forAll(v2Edges, v2EdgeI)
1369  {
1370  const edge& e = edges[v2Edges[v2EdgeI]];
1371 
1372  label vertI = e.otherVertex(v2);
1373 
1374  vertexNeighbours.insert(vertI);
1375  }
1376  return vertexNeighbours.toc();
1377 }
1378 
1379 
1381 //void Foam::triSurfaceTools::orderVertices
1382 //(
1383 // const labelledTri& f,
1384 // const label v1,
1385 // const label v2,
1386 // label& vA,
1387 // label& vB
1388 //)
1389 //{
1390 // // Order v1, v2 in anticlockwise order.
1391 // bool reverse = false;
1392 //
1393 // if (f[0] == v1)
1394 // {
1395 // if (f[1] != v2)
1396 // {
1397 // reverse = true;
1398 // }
1399 // }
1400 // else if (f[1] == v1)
1401 // {
1402 // if (f[2] != v2)
1403 // {
1404 // reverse = true;
1405 // }
1406 // }
1407 // else
1408 // {
1409 // if (f[0] != v2)
1410 // {
1411 // reverse = true;
1412 // }
1413 // }
1414 //
1415 // if (reverse)
1416 // {
1417 // vA = v2;
1418 // vB = v1;
1419 // }
1420 // else
1421 // {
1422 // vA = v1;
1423 // vB = v2;
1424 // }
1425 //}
1426 
1427 
1428 // Get the other face using edgeI
1431  const triSurface& surf,
1432  const label faceI,
1433  const label edgeI
1434 )
1435 {
1436  const labelList& myFaces = surf.edgeFaces()[edgeI];
1437 
1438  if (myFaces.size() != 2)
1439  {
1440  return -1;
1441  }
1442  else
1443  {
1444  if (faceI == myFaces[0])
1445  {
1446  return myFaces[1];
1447  }
1448  else
1449  {
1450  return myFaces[0];
1451  }
1452  }
1453 }
1454 
1455 
1456 // Get the two edges on faceI counterclockwise after edgeI
1459  const triSurface& surf,
1460  const label faceI,
1461  const label edgeI,
1462  label& e1,
1463  label& e2
1464 )
1465 {
1466  const labelList& eFaces = surf.faceEdges()[faceI];
1467 
1468  label i0 = findIndex(eFaces, edgeI);
1469 
1470  if (i0 == -1)
1471  {
1473  << "Edge " << surf.edges()[edgeI] << " not in face "
1474  << surf.localFaces()[faceI] << abort(FatalError);
1475  }
1476 
1477  label i1 = eFaces.fcIndex(i0);
1478  label i2 = eFaces.fcIndex(i1);
1479 
1480  e1 = eFaces[i1];
1481  e2 = eFaces[i2];
1482 }
1483 
1484 
1485 // Get the two vertices on faceI counterclockwise vertI
1488  const triSurface& surf,
1489  const label faceI,
1490  const label vertI,
1491  label& vert1I,
1492  label& vert2I
1493 )
1494 {
1495  const labelledTri& f = surf.localFaces()[faceI];
1496 
1497  if (vertI == f[0])
1498  {
1499  vert1I = f[1];
1500  vert2I = f[2];
1501  }
1502  else if (vertI == f[1])
1503  {
1504  vert1I = f[2];
1505  vert2I = f[0];
1506  }
1507  else if (vertI == f[2])
1508  {
1509  vert1I = f[0];
1510  vert2I = f[1];
1511  }
1512  else
1513  {
1515  << "Vertex " << vertI << " not in face " << f << abort(FatalError);
1516  }
1517 }
1518 
1519 
1520 // Get edge opposite vertex
1523  const triSurface& surf,
1524  const label faceI,
1525  const label vertI
1526 )
1527 {
1528  const labelList& myEdges = surf.faceEdges()[faceI];
1529 
1530  forAll(myEdges, myEdgeI)
1531  {
1532  label edgeI = myEdges[myEdgeI];
1533 
1534  const edge& e = surf.edges()[edgeI];
1535 
1536  if ((e.start() != vertI) && (e.end() != vertI))
1537  {
1538  return edgeI;
1539  }
1540  }
1541 
1543  << "Cannot find vertex " << vertI << " in edges of face " << faceI
1544  << abort(FatalError);
1545 
1546  return -1;
1547 }
1548 
1549 
1550 // Get vertex opposite edge
1553  const triSurface& surf,
1554  const label faceI,
1555  const label edgeI
1556 )
1557 {
1558  const triSurface::FaceType& f = surf.localFaces()[faceI];
1559  const edge& e = surf.edges()[edgeI];
1560 
1561  forAll(f, fp)
1562  {
1563  label vertI = f[fp];
1564 
1565  if (vertI != e.start() && vertI != e.end())
1566  {
1567  return vertI;
1568  }
1569  }
1570 
1572  << "Cannot find vertex opposite edge " << edgeI << " vertices " << e
1573  << " in face " << faceI << " vertices " << f << abort(FatalError);
1574 
1575  return -1;
1576 }
1577 
1578 
1579 // Returns edge label connecting v1, v2
1582  const triSurface& surf,
1583  const label v1,
1584  const label v2
1585 )
1586 {
1587  const labelList& v1Edges = surf.pointEdges()[v1];
1588 
1589  forAll(v1Edges, v1EdgeI)
1590  {
1591  label edgeI = v1Edges[v1EdgeI];
1592  const edge& e = surf.edges()[edgeI];
1593 
1594  if ((e.start() == v2) || (e.end() == v2))
1595  {
1596  return edgeI;
1597  }
1598  }
1599  return -1;
1600 }
1601 
1602 
1603 // Return index of triangle (or -1) using all three edges
1606  const triSurface& surf,
1607  const label e0I,
1608  const label e1I,
1609  const label e2I
1610 )
1611 {
1612  if ((e0I == e1I) || (e0I == e2I) || (e1I == e2I))
1613  {
1615  << "Duplicate edge labels : e0:" << e0I << " e1:" << e1I
1616  << " e2:" << e2I
1617  << abort(FatalError);
1618  }
1619 
1620  const labelList& eFaces = surf.edgeFaces()[e0I];
1621 
1622  forAll(eFaces, eFaceI)
1623  {
1624  label faceI = eFaces[eFaceI];
1625 
1626  const labelList& myEdges = surf.faceEdges()[faceI];
1627 
1628  if
1629  (
1630  (myEdges[0] == e1I)
1631  || (myEdges[1] == e1I)
1632  || (myEdges[2] == e1I)
1633  )
1634  {
1635  if
1636  (
1637  (myEdges[0] == e2I)
1638  || (myEdges[1] == e2I)
1639  || (myEdges[2] == e2I)
1640  )
1641  {
1642  return faceI;
1643  }
1644  }
1645  }
1646  return -1;
1647 }
1648 
1649 
1650 // Collapse indicated edges. Return new tri surface.
1653  const triSurface& surf,
1654  const labelList& collapsableEdges
1655 )
1656 {
1657  pointField edgeMids(surf.nEdges());
1658 
1659  forAll(edgeMids, edgeI)
1660  {
1661  const edge& e = surf.edges()[edgeI];
1662 
1663  edgeMids[edgeI] =
1664  0.5
1665  * (
1666  surf.localPoints()[e.start()]
1667  + surf.localPoints()[e.end()]
1668  );
1669  }
1670 
1671 
1672  labelList faceStatus(surf.size(), ANYEDGE);
1673 
1675  //forAll(edges, edgeI)
1676  //{
1677  // const labelList& neighbours = edgeFaces[edgeI];
1678  //
1679  // if ((neighbours.size() != 2) && (neighbours.size() != 1))
1680  // {
1681  // FatalErrorInFunction
1682  // << abort(FatalError);
1683  // }
1684  //
1685  // if (neighbours.size() == 2)
1686  // {
1687  // if (surf[neighbours[0]].region() != surf[neighbours[1]].region())
1688  // {
1689  // // Neighbours on different regions. For now dont allow
1690  // // any collapse.
1691  // //Pout<< "protecting face " << neighbours[0]
1692  // // << ' ' << neighbours[1] << endl;
1693  // faceStatus[neighbours[0]] = NOEDGE;
1694  // faceStatus[neighbours[1]] = NOEDGE;
1695  // }
1696  // }
1697  //}
1698 
1699  return collapseEdges(surf, collapsableEdges, edgeMids, faceStatus);
1700 }
1701 
1702 
1703 // Collapse indicated edges. Return new tri surface.
1706  const triSurface& surf,
1707  const labelList& collapseEdgeLabels,
1708  const pointField& edgeMids,
1709  labelList& faceStatus
1710 )
1711 {
1712  const labelListList& edgeFaces = surf.edgeFaces();
1713  const pointField& localPoints = surf.localPoints();
1714  const edgeList& edges = surf.edges();
1715 
1716  // Storage for new points
1717  pointField newPoints(localPoints);
1718 
1719  // Map for old to new points
1720  labelList pointMap(localPoints.size());
1721  forAll(localPoints, pointI)
1722  {
1723  pointMap[pointI] = pointI;
1724  }
1725 
1726 
1727  // Do actual 'collapsing' of edges
1728 
1729  forAll(collapseEdgeLabels, collapseEdgeI)
1730  {
1731  const label edgeI = collapseEdgeLabels[collapseEdgeI];
1732 
1733  if ((edgeI < 0) || (edgeI >= surf.nEdges()))
1734  {
1736  << "Edge label outside valid range." << endl
1737  << "edge label:" << edgeI << endl
1738  << "total number of edges:" << surf.nEdges() << endl
1739  << abort(FatalError);
1740  }
1741 
1742  const labelList& neighbours = edgeFaces[edgeI];
1743 
1744  if (neighbours.size() == 2)
1745  {
1746  const label stat0 = faceStatus[neighbours[0]];
1747  const label stat1 = faceStatus[neighbours[1]];
1748 
1749  // Check faceStatus to make sure this one can be collapsed
1750  if
1751  (
1752  ((stat0 == ANYEDGE) || (stat0 == edgeI))
1753  && ((stat1 == ANYEDGE) || (stat1 == edgeI))
1754  )
1755  {
1756  const edge& e = edges[edgeI];
1757 
1758  // Set up mapping to 'collapse' points of edge
1759  if
1760  (
1761  (pointMap[e.start()] != e.start())
1762  || (pointMap[e.end()] != e.end())
1763  )
1764  {
1766  << "points already mapped. Double collapse." << endl
1767  << "edgeI:" << edgeI
1768  << " start:" << e.start()
1769  << " end:" << e.end()
1770  << " pointMap[start]:" << pointMap[e.start()]
1771  << " pointMap[end]:" << pointMap[e.end()]
1772  << abort(FatalError);
1773  }
1774 
1775  const label minVert = min(e.start(), e.end());
1776  pointMap[e.start()] = minVert;
1777  pointMap[e.end()] = minVert;
1778 
1779  // Move shared vertex to mid of edge
1780  newPoints[minVert] = edgeMids[edgeI];
1781 
1782  // Protect neighbouring faces
1783  protectNeighbours(surf, e.start(), faceStatus);
1784  protectNeighbours(surf, e.end(), faceStatus);
1785  protectNeighbours
1786  (
1787  surf,
1788  oppositeVertex(surf, neighbours[0], edgeI),
1789  faceStatus
1790  );
1791  protectNeighbours
1792  (
1793  surf,
1794  oppositeVertex(surf, neighbours[1], edgeI),
1795  faceStatus
1796  );
1797 
1798  // Mark all collapsing faces
1799  labelList collapseFaces =
1800  getCollapsedFaces
1801  (
1802  surf,
1803  edgeI
1804  ).toc();
1805 
1806  forAll(collapseFaces, collapseI)
1807  {
1808  faceStatus[collapseFaces[collapseI]] = COLLAPSED;
1809  }
1810  }
1811  }
1812  }
1813 
1814 
1815  // Storage for new triangles
1816  List<labelledTri> newTris(surf.size());
1817  label newTriI = 0;
1818 
1819  const List<labelledTri>& localFaces = surf.localFaces();
1820 
1821 
1822  // Get only non-collapsed triangles and renumber vertex labels.
1823  forAll(localFaces, faceI)
1824  {
1825  const labelledTri& f = localFaces[faceI];
1826 
1827  const label a = pointMap[f[0]];
1828  const label b = pointMap[f[1]];
1829  const label c = pointMap[f[2]];
1830 
1831  if
1832  (
1833  (a != b) && (a != c) && (b != c)
1834  && (faceStatus[faceI] != COLLAPSED)
1835  )
1836  {
1837  // uncollapsed triangle
1838  newTris[newTriI++] = labelledTri(a, b, c, f.region());
1839  }
1840  else
1841  {
1842  //Pout<< "Collapsed triangle " << faceI
1843  // << " vertices:" << f << endl;
1844  }
1845  }
1846  newTris.setSize(newTriI);
1847 
1848 
1849 
1850  // Pack faces
1851 
1852  triSurface tempSurf(newTris, surf.patches(), newPoints);
1853 
1854  return
1855  triSurface
1856  (
1857  tempSurf.localFaces(),
1858  tempSurf.patches(),
1859  tempSurf.localPoints()
1860  );
1861 }
1862 
1863 
1864 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1865 
1868  const triSurface& surf,
1869  const labelList& refineFaces
1870 )
1871 {
1872  List<refineType> refineStatus(surf.size(), NONE);
1873 
1874  // Mark & propagate refinement
1875  forAll(refineFaces, refineFaceI)
1876  {
1877  calcRefineStatus(surf, refineFaces[refineFaceI], refineStatus);
1878  }
1879 
1880  // Do actual refinement
1881  return doRefine(surf, refineStatus);
1882 }
1883 
1884 
1887  const triSurface& surf,
1888  const labelList& refineEdges
1889 )
1890 {
1891  // Storage for marking faces
1892  List<refineType> refineStatus(surf.size(), NONE);
1893 
1894  // Storage for new faces
1895  DynamicList<labelledTri> newFaces(0);
1896 
1897  pointField newPoints(surf.localPoints());
1898  newPoints.setSize(surf.nPoints() + surf.nEdges());
1899  label newPointI = surf.nPoints();
1900 
1901 
1902  // Refine edges
1903  forAll(refineEdges, refineEdgeI)
1904  {
1905  label edgeI = refineEdges[refineEdgeI];
1906 
1907  const labelList& myFaces = surf.edgeFaces()[edgeI];
1908 
1909  bool neighbourIsRefined= false;
1910 
1911  forAll(myFaces, myFaceI)
1912  {
1913  if (refineStatus[myFaces[myFaceI]] != NONE)
1914  {
1915  neighbourIsRefined = true;
1916  }
1917  }
1918 
1919  // Only refine if none of the faces is refined
1920  if (!neighbourIsRefined)
1921  {
1922  // Refine edge
1923  const edge& e = surf.edges()[edgeI];
1924 
1925  point mid =
1926  0.5
1927  * (
1928  surf.localPoints()[e.start()]
1929  + surf.localPoints()[e.end()]
1930  );
1931 
1932  newPoints[newPointI] = mid;
1933 
1934  // Refine faces using edge
1935  forAll(myFaces, myFaceI)
1936  {
1937  // Add faces to newFaces
1938  greenRefine
1939  (
1940  surf,
1941  myFaces[myFaceI],
1942  edgeI,
1943  newPointI,
1944  newFaces
1945  );
1946 
1947  // Mark as refined
1948  refineStatus[myFaces[myFaceI]] = GREEN;
1949  }
1950 
1951  newPointI++;
1952  }
1953  }
1954 
1955  // Add unrefined faces
1956  forAll(surf.localFaces(), faceI)
1957  {
1958  if (refineStatus[faceI] == NONE)
1959  {
1960  newFaces.append(surf.localFaces()[faceI]);
1961  }
1962  }
1963 
1964  newFaces.shrink();
1965  newPoints.setSize(newPointI);
1966 
1967  return triSurface(newFaces, surf.patches(), newPoints, true);
1968 }
1969 
1970 
1971 
1972 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1973 // Geometric helper functions:
1974 
1975 
1976 // Returns element in edgeIndices with minimum length
1979  const triSurface& surf,
1980  const labelList& edgeIndices
1981 )
1982 {
1983  scalar minLength = GREAT;
1984  label minIndex = -1;
1985  forAll(edgeIndices, i)
1986  {
1987  const edge& e = surf.edges()[edgeIndices[i]];
1988 
1989  scalar length =
1990  mag
1991  (
1992  surf.localPoints()[e.end()]
1993  - surf.localPoints()[e.start()]
1994  );
1995 
1996  if (length < minLength)
1997  {
1998  minLength = length;
1999  minIndex = i;
2000  }
2001  }
2002  return edgeIndices[minIndex];
2003 }
2004 
2005 
2006 // Returns element in edgeIndices with maximum length
2009  const triSurface& surf,
2010  const labelList& edgeIndices
2011 )
2012 {
2013  scalar maxLength = -GREAT;
2014  label maxIndex = -1;
2015  forAll(edgeIndices, i)
2016  {
2017  const edge& e = surf.edges()[edgeIndices[i]];
2018 
2019  scalar length =
2020  mag
2021  (
2022  surf.localPoints()[e.end()]
2023  - surf.localPoints()[e.start()]
2024  );
2025 
2026  if (length > maxLength)
2027  {
2028  maxLength = length;
2029  maxIndex = i;
2030  }
2031  }
2032  return edgeIndices[maxIndex];
2033 }
2034 
2035 
2036 // Merge points and reconstruct surface
2039  const triSurface& surf,
2040  const scalar mergeTol
2041 )
2042 {
2043  pointField newPoints(surf.nPoints());
2044 
2045  labelList pointMap(surf.nPoints());
2046 
2047  bool hasMerged = Foam::mergePoints
2048  (
2049  surf.localPoints(),
2050  mergeTol,
2051  false,
2052  pointMap,
2053  newPoints
2054  );
2055 
2056  if (hasMerged)
2057  {
2058  // Pack the triangles
2059 
2060  // Storage for new triangles
2061  List<labelledTri> newTriangles(surf.size());
2062  label newTriangleI = 0;
2063 
2064  forAll(surf, faceI)
2065  {
2066  const labelledTri& f = surf.localFaces()[faceI];
2067 
2068  label newA = pointMap[f[0]];
2069  label newB = pointMap[f[1]];
2070  label newC = pointMap[f[2]];
2071 
2072  if ((newA != newB) && (newA != newC) && (newB != newC))
2073  {
2074  newTriangles[newTriangleI++] =
2075  labelledTri(newA, newB, newC, f.region());
2076  }
2077  }
2078  newTriangles.setSize(newTriangleI);
2079 
2080  return triSurface
2081  (
2082  newTriangles,
2083  surf.patches(),
2084  newPoints,
2085  true //reuse storage
2086  );
2087  }
2088  else
2089  {
2090  return surf;
2091  }
2092 }
2093 
2094 
2095 // Calculate normal on triangle
2098  const triSurface& surf,
2099  const label nearestFaceI,
2100  const point& nearestPt
2101 )
2102 {
2103  const triSurface::FaceType& f = surf[nearestFaceI];
2104  const pointField& points = surf.points();
2105 
2106  label nearType, nearLabel;
2107 
2108  f.nearestPointClassify(nearestPt, points, nearType, nearLabel);
2109 
2110  if (nearType == triPointRef::NONE)
2111  {
2112  // Nearest to face
2113  return surf.faceNormals()[nearestFaceI];
2114  }
2115  else if (nearType == triPointRef::EDGE)
2116  {
2117  // Nearest to edge. Assume order of faceEdges same as face vertices.
2118  label edgeI = surf.faceEdges()[nearestFaceI][nearLabel];
2119 
2120  // Calculate edge normal by averaging face normals
2121  const labelList& eFaces = surf.edgeFaces()[edgeI];
2122 
2123  vector edgeNormal(vector::zero);
2124 
2125  forAll(eFaces, i)
2126  {
2127  edgeNormal += surf.faceNormals()[eFaces[i]];
2128  }
2129  return edgeNormal/(mag(edgeNormal) + VSMALL);
2130  }
2131  else
2132  {
2133  // Nearest to point
2134  const triSurface::FaceType& localF = surf.localFaces()[nearestFaceI];
2135  return surf.pointNormals()[localF[nearLabel]];
2136  }
2137 }
2138 
2139 
2142  const triSurface& surf,
2143  const point& sample,
2144  const point& nearestPoint,
2145  const label edgeI
2146 )
2147 {
2148  const labelList& eFaces = surf.edgeFaces()[edgeI];
2149 
2150  if (eFaces.size() != 2)
2151  {
2152  // Surface not closed.
2153  return UNKNOWN;
2154  }
2155  else
2156  {
2157  const vectorField& faceNormals = surf.faceNormals();
2158 
2159  // Compare to bisector. This is actually correct since edge is
2160  // nearest so there is a knife-edge.
2161 
2162  vector n = 0.5*(faceNormals[eFaces[0]] + faceNormals[eFaces[1]]);
2163 
2164  if (((sample - nearestPoint) & n) > 0)
2165  {
2166  return OUTSIDE;
2167  }
2168  else
2169  {
2170  return INSIDE;
2171  }
2172  }
2173 }
2174 
2175 
2176 // Calculate normal on triangle
2179  const triSurface& surf,
2180  const point& sample,
2181  const label nearestFaceI
2182 )
2183 {
2184  const triSurface::FaceType& f = surf[nearestFaceI];
2185  const pointField& points = surf.points();
2186 
2187  // Find where point is on face
2188  label nearType, nearLabel;
2189 
2190  pointHit pHit = f.nearestPointClassify(sample, points, nearType, nearLabel);
2191 
2192  const point& nearestPoint(pHit.rawPoint());
2193 
2194  if (nearType == triPointRef::NONE)
2195  {
2196  vector sampleNearestVec = (sample - nearestPoint);
2197 
2198  // Nearest to face interior. Use faceNormal to determine side
2199  scalar c = sampleNearestVec & surf.faceNormals()[nearestFaceI];
2200 
2201  // // If the sample is essentially on the face, do not check for
2202  // // it being perpendicular.
2203 
2204  // scalar magSampleNearestVec = mag(sampleNearestVec);
2205 
2206  // if (magSampleNearestVec > SMALL)
2207  // {
2208  // c /= magSampleNearestVec*mag(surf.faceNormals()[nearestFaceI]);
2209 
2210  // if (mag(c) < 0.99)
2211  // {
2212  // FatalErrorInFunction
2213  // << "nearestPoint identified as being on triangle face "
2214  // << "but vector from nearestPoint to sample is not "
2215  // << "perpendicular to the normal." << nl
2216  // << "sample: " << sample << nl
2217  // << "nearestPoint: " << nearestPoint << nl
2218  // << "sample - nearestPoint: "
2219  // << sample - nearestPoint << nl
2220  // << "normal: " << surf.faceNormals()[nearestFaceI] << nl
2221  // << "mag(sample - nearestPoint): "
2222  // << mag(sample - nearestPoint) << nl
2223  // << "normalised dot product: " << c << nl
2224  // << "triangle vertices: " << nl
2225  // << " " << points[f[0]] << nl
2226  // << " " << points[f[1]] << nl
2227  // << " " << points[f[2]] << nl
2228  // << abort(FatalError);
2229  // }
2230  // }
2231 
2232  if (c > 0)
2233  {
2234  return OUTSIDE;
2235  }
2236  else
2237  {
2238  return INSIDE;
2239  }
2240  }
2241  else if (nearType == triPointRef::EDGE)
2242  {
2243  // Nearest to edge nearLabel. Note that this can only be a knife-edge
2244  // situation since otherwise the nearest point could never be the edge.
2245 
2246  // Get the edge. Assume order of faceEdges same as face vertices.
2247  label edgeI = surf.faceEdges()[nearestFaceI][nearLabel];
2248 
2249  // if (debug)
2250  // {
2251  // // Check order of faceEdges same as face vertices.
2252  // const edge& e = surf.edges()[edgeI];
2253  // const labelList& meshPoints = surf.meshPoints();
2254  // const edge meshEdge(meshPoints[e[0]], meshPoints[e[1]]);
2255 
2256  // if
2257  // (
2258  // meshEdge
2259  // != edge(f[nearLabel], f[f.fcIndex(nearLabel)])
2260  // )
2261  // {
2262  // FatalErrorInFunction
2263  // << "Edge:" << edgeI << " local vertices:" << e
2264  // << " mesh vertices:" << meshEdge
2265  // << " not at position " << nearLabel
2266  // << " in face " << f
2267  // << abort(FatalError);
2268  // }
2269  // }
2270 
2271  return edgeSide(surf, sample, nearestPoint, edgeI);
2272  }
2273  else
2274  {
2275  // Nearest to point. Could use pointNormal here but is not correct.
2276  // Instead determine which edge using point is nearest and use test
2277  // above (nearType == triPointRef::EDGE).
2278 
2279 
2280  const triSurface::FaceType& localF = surf.localFaces()[nearestFaceI];
2281  label nearPointI = localF[nearLabel];
2282 
2283  const edgeList& edges = surf.edges();
2284  const pointField& localPoints = surf.localPoints();
2285  const point& base = localPoints[nearPointI];
2286 
2287  const labelList& pEdges = surf.pointEdges()[nearPointI];
2288 
2289  scalar minDistSqr = Foam::sqr(GREAT);
2290  label minEdgeI = -1;
2291 
2292  forAll(pEdges, i)
2293  {
2294  label edgeI = pEdges[i];
2295 
2296  const edge& e = edges[edgeI];
2297 
2298  label otherPointI = e.otherVertex(nearPointI);
2299 
2300  // Get edge normal.
2301  vector eVec(localPoints[otherPointI] - base);
2302  scalar magEVec = mag(eVec);
2303 
2304  if (magEVec > VSMALL)
2305  {
2306  eVec /= magEVec;
2307 
2308  // Get point along vector and determine closest.
2309  const point perturbPoint = base + eVec;
2310 
2311  scalar distSqr = Foam::magSqr(sample - perturbPoint);
2312 
2313  if (distSqr < minDistSqr)
2314  {
2315  minDistSqr = distSqr;
2316  minEdgeI = edgeI;
2317  }
2318  }
2319  }
2320 
2321  if (minEdgeI == -1)
2322  {
2324  << "Problem: did not find edge closer than " << minDistSqr
2325  << abort(FatalError);
2326  }
2327 
2328  return edgeSide(surf, sample, nearestPoint, minEdgeI);
2329  }
2330 }
2331 
2332 
2333 // triangulation of boundaryMesh
2336  const polyBoundaryMesh& bMesh,
2337  const labelHashSet& includePatches,
2338  const bool verbose
2339 )
2340 {
2341  const polyMesh& mesh = bMesh.mesh();
2342 
2343  // Storage for surfaceMesh. Size estimate.
2344  DynamicList<labelledTri> triangles
2345  (
2346  mesh.nFaces() - mesh.nInternalFaces()
2347  );
2348 
2349  label newPatchI = 0;
2350 
2351  forAllConstIter(labelHashSet, includePatches, iter)
2352  {
2353  const label patchI = iter.key();
2354  const polyPatch& patch = bMesh[patchI];
2355  const pointField& points = patch.points();
2356 
2357  label nTriTotal = 0;
2358 
2359  forAll(patch, patchFaceI)
2360  {
2361  const face& f = patch[patchFaceI];
2362 
2363  faceList triFaces(f.nTriangles(points));
2364 
2365  label nTri = 0;
2366 
2367  f.triangles(points, nTri, triFaces);
2368 
2369  forAll(triFaces, triFaceI)
2370  {
2371  const face& f = triFaces[triFaceI];
2372 
2373  triangles.append(labelledTri(f[0], f[1], f[2], newPatchI));
2374 
2375  nTriTotal++;
2376  }
2377  }
2378 
2379  if (verbose)
2380  {
2381  Pout<< patch.name() << " : generated " << nTriTotal
2382  << " triangles from " << patch.size() << " faces with"
2383  << " new patchid " << newPatchI << endl;
2384  }
2385 
2386  newPatchI++;
2387  }
2388  triangles.shrink();
2389 
2390  // Create globally numbered tri surface
2391  triSurface rawSurface(triangles, mesh.points());
2392 
2393  // Create locally numbered tri surface
2395  (
2396  rawSurface.localFaces(),
2397  rawSurface.localPoints()
2398  );
2399 
2400  // Add patch names to surface
2401  surface.patches().setSize(newPatchI);
2402 
2403  newPatchI = 0;
2404 
2405  forAllConstIter(labelHashSet, includePatches, iter)
2406  {
2407  const label patchI = iter.key();
2408  const polyPatch& patch = bMesh[patchI];
2409 
2410  surface.patches()[newPatchI].name() = patch.name();
2411  surface.patches()[newPatchI].geometricType() = patch.type();
2412 
2413  newPatchI++;
2414  }
2415 
2416  return surface;
2417 }
2418 
2419 
2420 // triangulation of boundaryMesh
2423  const polyBoundaryMesh& bMesh,
2424  const labelHashSet& includePatches,
2425  const bool verbose
2426 )
2427 {
2428  const polyMesh& mesh = bMesh.mesh();
2429 
2430  // Storage for new points = meshpoints + face centres.
2431  const pointField& points = mesh.points();
2432  const pointField& faceCentres = mesh.faceCentres();
2433 
2434  pointField newPoints(points.size() + faceCentres.size());
2435 
2436  label newPointI = 0;
2437 
2438  forAll(points, pointI)
2439  {
2440  newPoints[newPointI++] = points[pointI];
2441  }
2442  forAll(faceCentres, faceI)
2443  {
2444  newPoints[newPointI++] = faceCentres[faceI];
2445  }
2446 
2447 
2448  // Count number of faces.
2449  DynamicList<labelledTri> triangles
2450  (
2451  mesh.nFaces() - mesh.nInternalFaces()
2452  );
2453 
2454  label newPatchI = 0;
2455 
2456  forAllConstIter(labelHashSet, includePatches, iter)
2457  {
2458  const label patchI = iter.key();
2459  const polyPatch& patch = bMesh[patchI];
2460 
2461  label nTriTotal = 0;
2462 
2463  forAll(patch, patchFaceI)
2464  {
2465  // Face in global coords.
2466  const face& f = patch[patchFaceI];
2467 
2468  // Index in newPointI of face centre.
2469  label fc = points.size() + patchFaceI + patch.start();
2470 
2471  forAll(f, fp)
2472  {
2473  label fp1 = f.fcIndex(fp);
2474 
2475  triangles.append(labelledTri(f[fp], f[fp1], fc, newPatchI));
2476 
2477  nTriTotal++;
2478  }
2479  }
2480 
2481  if (verbose)
2482  {
2483  Pout<< patch.name() << " : generated " << nTriTotal
2484  << " triangles from " << patch.size() << " faces with"
2485  << " new patchid " << newPatchI << endl;
2486  }
2487 
2488  newPatchI++;
2489  }
2490  triangles.shrink();
2491 
2492 
2493  // Create globally numbered tri surface
2494  triSurface rawSurface(triangles, newPoints);
2495 
2496  // Create locally numbered tri surface
2498  (
2499  rawSurface.localFaces(),
2500  rawSurface.localPoints()
2501  );
2502 
2503  // Add patch names to surface
2504  surface.patches().setSize(newPatchI);
2505 
2506  newPatchI = 0;
2507 
2508  forAllConstIter(labelHashSet, includePatches, iter)
2509  {
2510  const label patchI = iter.key();
2511  const polyPatch& patch = bMesh[patchI];
2512 
2513  surface.patches()[newPatchI].name() = patch.name();
2514  surface.patches()[newPatchI].geometricType() = patch.type();
2515 
2516  newPatchI++;
2517  }
2518 
2519  return surface;
2520 }
2521 
2522 
2524 {
2525  // Vertices in geompack notation. Note that could probably just use
2526  // pts.begin() if double precision.
2527  List<doubleScalar> geompackVertices(2*pts.size());
2528  label doubleI = 0;
2529  forAll(pts, i)
2530  {
2531  geompackVertices[doubleI++] = pts[i][0];
2532  geompackVertices[doubleI++] = pts[i][1];
2533  }
2534 
2535  // Storage for triangles
2536  int m2 = 3;
2537  List<int> triangle_node(m2*3*pts.size());
2538  List<int> triangle_neighbor(m2*3*pts.size());
2539 
2540  // Triangulate
2541  int nTris = 0;
2542  int err = dtris2
2543  (
2544  pts.size(),
2545  geompackVertices.begin(),
2546  &nTris,
2547  triangle_node.begin(),
2548  triangle_neighbor.begin()
2549  );
2550 
2551  if (err != 0)
2552  {
2554  << "Failed dtris2 with vertices:" << pts.size()
2555  << abort(FatalError);
2556  }
2557 
2558  // Trim
2559  triangle_node.setSize(3*nTris);
2560  triangle_neighbor.setSize(3*nTris);
2561 
2562  // Convert to triSurface.
2563  List<labelledTri> faces(nTris);
2564 
2565  forAll(faces, i)
2566  {
2567  faces[i] = labelledTri
2568  (
2569  triangle_node[3*i]-1,
2570  triangle_node[3*i+1]-1,
2571  triangle_node[3*i+2]-1,
2572  0
2573  );
2574  }
2575 
2576  pointField points(pts.size());
2577  forAll(pts, i)
2578  {
2579  points[i][0] = pts[i][0];
2580  points[i][1] = pts[i][1];
2581  points[i][2] = 0.0;
2582  }
2583 
2584  return triSurface(faces, points);
2585 }
2586 
2587 
2590  const triPointRef& tri,
2591  const point& p,
2592  FixedList<scalar, 3>& weights
2593 )
2594 {
2595  // calculate triangle edge vectors and triangle face normal
2596  // the 'i':th edge is opposite node i
2598  edge[0] = tri.c()-tri.b();
2599  edge[1] = tri.a()-tri.c();
2600  edge[2] = tri.b()-tri.a();
2601 
2602  vector triangleFaceNormal = edge[1] ^ edge[2];
2603 
2604  // calculate edge normal (pointing inwards)
2606  for (label i=0; i<3; i++)
2607  {
2608  normal[i] = triangleFaceNormal ^ edge[i];
2609  normal[i] /= mag(normal[i]) + VSMALL;
2610  }
2611 
2612  weights[0] = ((p-tri.b()) & normal[0]) / max(VSMALL, normal[0] & edge[1]);
2613  weights[1] = ((p-tri.c()) & normal[1]) / max(VSMALL, normal[1] & edge[2]);
2614  weights[2] = ((p-tri.a()) & normal[2]) / max(VSMALL, normal[2] & edge[0]);
2615 }
2616 
2617 
2618 // Calculate weighting factors from samplePts to triangle it is in.
2619 // Uses linear search.
2622  const triSurface& s,
2623  const pointField& samplePts,
2624  List<FixedList<label, 3> >& allVerts,
2625  List<FixedList<scalar, 3> >& allWeights
2626 )
2627 {
2628  allVerts.setSize(samplePts.size());
2629  allWeights.setSize(samplePts.size());
2630 
2631  const pointField& points = s.points();
2632 
2633  forAll(samplePts, i)
2634  {
2635  const point& samplePt = samplePts[i];
2636 
2637 
2638  FixedList<label, 3>& verts = allVerts[i];
2639  FixedList<scalar, 3>& weights = allWeights[i];
2640 
2641  scalar minDistance = GREAT;
2642 
2643  forAll(s, faceI)
2644  {
2645  const labelledTri& f = s[faceI];
2646 
2647  triPointRef tri(f.tri(points));
2648 
2649  label nearType, nearLabel;
2650 
2651  pointHit nearest = tri.nearestPointClassify
2652  (
2653  samplePt,
2654  nearType,
2655  nearLabel
2656  );
2657 
2658  if (nearest.hit())
2659  {
2660  // samplePt inside triangle
2661  verts[0] = f[0];
2662  verts[1] = f[1];
2663  verts[2] = f[2];
2664 
2665  calcInterpolationWeights(tri, nearest.rawPoint(), weights);
2666 
2667  //Pout<< "calcScalingFactors : samplePt:" << samplePt
2668  // << " inside triangle:" << faceI
2669  // << " verts:" << verts
2670  // << " weights:" << weights
2671  // << endl;
2672 
2673  break;
2674  }
2675  else if (nearest.distance() < minDistance)
2676  {
2677  minDistance = nearest.distance();
2678 
2679  // Outside triangle. Store nearest.
2680 
2681  if (nearType == triPointRef::POINT)
2682  {
2683  verts[0] = f[nearLabel];
2684  weights[0] = 1;
2685  verts[1] = -1;
2686  weights[1] = -GREAT;
2687  verts[2] = -1;
2688  weights[2] = -GREAT;
2689 
2690  //Pout<< "calcScalingFactors : samplePt:" << samplePt
2691  // << " distance:" << nearest.distance()
2692  // << " from point:" << points[f[nearLabel]]
2693  // << endl;
2694  }
2695  else if (nearType == triPointRef::EDGE)
2696  {
2697  verts[0] = f[nearLabel];
2698  verts[1] = f[f.fcIndex(nearLabel)];
2699  verts[2] = -1;
2700 
2701  const point& p0 = points[verts[0]];
2702  const point& p1 = points[verts[1]];
2703 
2704  scalar s = min
2705  (
2706  1,
2707  max
2708  (
2709  0,
2710  mag(nearest.rawPoint() - p0)/mag(p1 - p0)
2711  )
2712  );
2713 
2714  // Interpolate
2715  weights[0] = 1 - s;
2716  weights[1] = s;
2717  weights[2] = -GREAT;
2718 
2719  //Pout<< "calcScalingFactors : samplePt:" << samplePt
2720  // << " distance:" << nearest.distance()
2721  // << " from edge:" << p0 << p1 << " s:" << s
2722  // << endl;
2723  }
2724  else
2725  {
2726  // triangle. Can only happen because of truncation errors.
2727  verts[0] = f[0];
2728  verts[1] = f[1];
2729  verts[2] = f[2];
2730 
2731  calcInterpolationWeights(tri, nearest.rawPoint(), weights);
2732 
2733  //Pout<< "calcScalingFactors : samplePt:" << samplePt
2734  // << " distance:" << nearest.distance()
2735  // << " to verts:" << verts
2736  // << " weights:" << weights
2737  // << endl;
2738 
2739  break;
2740  }
2741  }
2742  }
2743  }
2744 }
2745 
2746 
2747 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2748 // Tracking:
2749 
2750 
2751 // Test point on surface to see if is on face,edge or point.
2754  const triSurface& s,
2755  const label triI,
2756  const point& trianglePoint
2757 )
2758 {
2759  surfaceLocation nearest;
2760 
2761  // Nearest point could be on point or edge. Retest.
2762  label index, elemType;
2763  //bool inside =
2764  triPointRef(s[triI].tri(s.points())).classify
2765  (
2766  trianglePoint,
2767  elemType,
2768  index
2769  );
2770 
2771  nearest.setPoint(trianglePoint);
2772 
2773  if (elemType == triPointRef::NONE)
2774  {
2775  nearest.setHit();
2776  nearest.setIndex(triI);
2777  nearest.elementType() = triPointRef::NONE;
2778  }
2779  else if (elemType == triPointRef::EDGE)
2780  {
2781  nearest.setMiss();
2782  nearest.setIndex(s.faceEdges()[triI][index]);
2783  nearest.elementType() = triPointRef::EDGE;
2784  }
2785  else //if (elemType == triPointRef::POINT)
2786  {
2787  nearest.setMiss();
2788  nearest.setIndex(s.localFaces()[triI][index]);
2789  nearest.elementType() = triPointRef::POINT;
2790  }
2791 
2792  return nearest;
2793 }
2794 
2795 
2798  const triSurface& s,
2799  const surfaceLocation& start,
2800  const surfaceLocation& end,
2801  const plane& cutPlane
2802 )
2803 {
2804  // Start off from starting point
2805  surfaceLocation nearest = start;
2806  nearest.setMiss();
2807 
2808  // See if in same triangle as endpoint. If so snap.
2809  snapToEnd(s, end, nearest);
2810 
2811  if (!nearest.hit())
2812  {
2813  // Not yet at end point
2814 
2815  if (start.elementType() == triPointRef::NONE)
2816  {
2817  // Start point is inside triangle. Trivial cases already handled
2818  // above.
2819 
2820  // end point is on edge or point so cross currrent triangle to
2821  // see which edge is cut.
2822 
2823  nearest = cutEdge
2824  (
2825  s,
2826  start.index(), // triangle
2827  -1, // excludeEdge
2828  -1, // excludePoint
2829  start.rawPoint(),
2830  cutPlane,
2831  end.rawPoint()
2832  );
2833  nearest.elementType() = triPointRef::EDGE;
2834  nearest.triangle() = start.index();
2835  nearest.setMiss();
2836  }
2837  else if (start.elementType() == triPointRef::EDGE)
2838  {
2839  // Pick connected triangle that is most in direction.
2840  const labelList& eFaces = s.edgeFaces()[start.index()];
2841 
2842  nearest = visitFaces
2843  (
2844  s,
2845  eFaces,
2846  start,
2847  start.index(), // excludeEdgeI
2848  -1, // excludePointI
2849  end,
2850  cutPlane
2851  );
2852  }
2853  else // start.elementType() == triPointRef::POINT
2854  {
2855  const labelList& pFaces = s.pointFaces()[start.index()];
2856 
2857  nearest = visitFaces
2858  (
2859  s,
2860  pFaces,
2861  start,
2862  -1, // excludeEdgeI
2863  start.index(), // excludePointI
2864  end,
2865  cutPlane
2866  );
2867  }
2868  snapToEnd(s, end, nearest);
2869  }
2870  return nearest;
2871 }
2872 
2873 
2876  const triSurface& s,
2877  const surfaceLocation& endInfo,
2878  const plane& cutPlane,
2879  surfaceLocation& hitInfo
2880 )
2881 {
2882  //OFstream str("track.obj");
2883  //label vertI = 0;
2884  //meshTools::writeOBJ(str, hitInfo.rawPoint());
2885  //vertI++;
2886 
2887  // Track across surface.
2888  while (true)
2889  {
2890  //Pout<< "Tracking from:" << nl
2891  // << " " << hitInfo.info()
2892  // << endl;
2893 
2894  hitInfo = trackToEdge
2895  (
2896  s,
2897  hitInfo,
2898  endInfo,
2899  cutPlane
2900  );
2901 
2902  //meshTools::writeOBJ(str, hitInfo.rawPoint());
2903  //vertI++;
2904  //str<< "l " << vertI-1 << ' ' << vertI << nl;
2905 
2906  //Pout<< "Tracked to:" << nl
2907  // << " " << hitInfo.info() << endl;
2908 
2909  if (hitInfo.hit() || hitInfo.triangle() == -1)
2910  {
2911  break;
2912  }
2913  }
2914 }
2915 
2916 
2917 // ************************************************************************* //
Foam::triSurfaceTools::triangulate
static triSurface triangulate(const polyBoundaryMesh &mBesh, const labelHashSet &includePatches, const bool verbose=false)
Simple triangulation of (selected patches of) boundaryMesh. Needs.
Definition: triSurfaceTools.C:2335
Foam::triSurfaceTools::writeOBJ
static void writeOBJ(const fileName &fName, const pointField &pts)
Write pointField to OBJ format file.
Definition: triSurfaceTools.C:1253
Foam::plane::normal
const vector & normal() const
Return plane normal.
Definition: plane.C:248
Foam::triSurfaceTools::collapseCreatesFold
bool collapseCreatesFold(const triSurface &surf, const label v1, const point &pt, const labelHashSet &collapsedFaces, const HashTable< label, label, Hash< label > > &edgeToEdge, const HashTable< label, label, Hash< label > > &edgeToFace, const scalar minCos)
Like collapseMinCosAngle but return true for value < minCos.
Definition: triSurfaceTools.C:684
Foam::PrimitivePatch::pointFaces
const labelListList & pointFaces() const
Return point-face addressing.
Definition: PrimitivePatchTemplate.C:352
Foam::PointIndexHit::index
label index() const
Return index.
Definition: PointIndexHit.H:111
Foam::triSurfaceTools::getEdge
static label getEdge(const triSurface &surf, const label vert1I, const label vert2I)
Returns edge label connecting v1, v2 (local numbering)
Definition: triSurfaceTools.C:1581
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Return edge-face addressing.
Definition: PrimitivePatchTemplate.C:292
Foam::PointHit::hit
bool hit() const
Is there a hit.
Definition: PointHit.H:120
Foam::triSurfaceTools::getMergedEdges
static void getMergedEdges(const triSurface &surf, const label edgeI, const labelHashSet &collapsedFaces, HashTable< label, label, Hash< label > > &edgeToEdge, HashTable< label, label, Hash< label > > &edgeToFace)
Definition: triSurfaceTools.C:450
p
p
Definition: pEqn.H:62
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::PointIndexHit::setIndex
void setIndex(const label index)
Definition: PointIndexHit.H:168
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatchTemplate.C:212
Foam::HashTable::toc
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:201
Foam::triSurfaceTools::otherFace
static label otherFace(const triSurface &surf, const label faceI, const label edgeI)
Get face connected to edge not faceI.
Definition: triSurfaceTools.C:1430
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::PointHit
This class describes the interaction of a face and a point. It carries the info of a successful hit a...
Definition: PointHit.H:51
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
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::triSurfaceTools::triangulateFaceCentre
triSurface triangulateFaceCentre(const polyBoundaryMesh &mBesh, const labelHashSet &includePatches, const bool verbose=false)
Face-centre triangulation of (selected patches of) boundaryMesh.
Definition: triSurfaceTools.C:2422
Foam::triSurfaceTools::maxEdge
static label maxEdge(const triSurface &surf, const labelList &edgeIndices)
Returns element in edgeIndices with minimum length.
Definition: triSurfaceTools.C:2008
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
Foam::triSurfaceTools::delaunay2D
static triSurface delaunay2D(const List< vector2D > &)
Do unconstrained Delaunay of points. Returns triSurface with 3D.
Definition: triSurfaceTools.C:2523
Foam::triSurfaceTools::ANYEDGE
static const label ANYEDGE
Face collapse status.
Definition: triSurfaceTools.H:366
Foam::triSurfaceTools::collapseEdges
static triSurface collapseEdges(const triSurface &surf, const labelList &collapsableEdges)
Create new triSurface by collapsing edges to edge mids.
Definition: triSurfaceTools.C:1652
Foam::PointHit::rawPoint
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
Foam::PrimitivePatch::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatchTemplate.H:299
Foam::triSurfaceTools::faceCosAngle
static scalar faceCosAngle(const point &pStart, const point &pEnd, const point &pLeft, const point &pRight)
Definition: triSurfaceTools.C:307
Foam::triSurfaceTools::mergePoints
static triSurface mergePoints(const triSurface &surf, const scalar mergeTol)
Merge points within distance.
Definition: triSurfaceTools.C:2038
Foam::surfaceLocation
Contains information about location on a triSurface:
Definition: surfaceLocation.H:60
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::PointIndexHit::hit
bool hit() const
Is there a hit.
Definition: PointIndexHit.H:105
dtris2
int dtris2(int point_num, double point_xy[], int *tri_num, int tri_vert[], int tri_nabe[])
Definition: geompack.C:949
Foam::triSurfaceTools::otherVertices
static void otherVertices(const triSurface &surf, const label faceI, const label vertI, label &vert1I, label &vert2I)
Get the two vertices (local numbering) on faceI counterclockwise.
Definition: triSurfaceTools.C:1487
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::PointIndexHit::setMiss
void setMiss()
Definition: PointIndexHit.H:158
Foam::PrimitivePatch::faceEdges
const labelListList & faceEdges() const
Return face-edge addressing.
Definition: PrimitivePatchTemplate.C:312
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::triSurfaceTools::doRefine
static triSurface doRefine(const triSurface &surf, const List< refineType > &refineStatus)
Definition: triSurfaceTools.C:160
polyMesh.H
Foam::triSurface::patches
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:301
Foam::HashSet< label, Hash< label > >
Foam::surfaceLocation::elementType
triPointRef::proxType & elementType()
Definition: surfaceLocation.H:107
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A
simpleMatrix< scalar > A(Nc)
Foam::triSurfaceTools::snapToEnd
static void snapToEnd(const triSurface &s, const surfaceLocation &endInfo, surfaceLocation &current)
Checks if current is on the same triangle as the endpoint.
Definition: triSurfaceTools.C:1036
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
OFstream.H
Foam::plane
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
Foam::triSurfaceTools::COLLAPSED
static const label COLLAPSED
Definition: triSurfaceTools.H:368
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:59
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::triSurfaceTools::calcRefineStatus
static void calcRefineStatus(const triSurface &surf, const label faceI, List< refineType > &refine)
Definition: triSurfaceTools.C:54
Foam::PointIndexHit::rawPoint
const Point & rawPoint() const
Return point with no checking.
Definition: PointIndexHit.H:143
Foam::triSurfaceTools::getTriangle
static label getTriangle(const triSurface &surf, const label e0I, const label e1I, const label e2I)
Return index of triangle (or -1) using all three edges.
Definition: triSurfaceTools.C:1605
Foam::Hash
Hash function class for primitives. All non-primitives used to hash entries on hash tables likely nee...
Definition: Hash.H:56
geompack.H
Foam::PointHit::distance
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
Foam::triangle::c
const Point & c() const
Return third vertex.
Definition: triangleI.H:95
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::triSurfaceTools::protectNeighbours
static void protectNeighbours(const triSurface &surf, const label vertI, labelList &faceStatus)
Definition: triSurfaceTools.C:333
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
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:2097
plane.H
Foam::triSurfaceTools::oppositeEdge
static label oppositeEdge(const triSurface &surf, const label faceI, const label vertI)
Get edge opposite vertex (local numbering)
Definition: triSurfaceTools.C:1522
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::triangle::a
const Point & a() const
Return first vertex.
Definition: triangleI.H:83
Foam::PrimitivePatch::pointEdges
const labelListList & pointEdges() const
Return point-edge addressing.
Definition: PrimitivePatchTemplate.C:332
Foam::DynamicList::clear
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:242
Foam::triSurfaceTools::redGreenRefine
static triSurface redGreenRefine(const triSurface &surf, const labelList &refineFaces)
Refine face by splitting all edges. Neighbouring face is.
Definition: triSurfaceTools.C:1867
Foam::triangle::NONE
@ NONE
Definition: triangle.H:96
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:2178
Foam::PointIndexHit::setHit
void setHit()
Definition: PointIndexHit.H:153
Foam::DynamicList::shrink
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:258
Foam::triSurfaceTools::oppositeVertex
static label oppositeVertex(const triSurface &surf, const label faceI, const label edgeI)
Get vertex (local numbering) opposite edge.
Definition: triSurfaceTools.C:1552
Foam::triSurfaceTools::edgeCosAngle
static scalar edgeCosAngle(const triSurface &surf, const label v1, const point &pt, const labelHashSet &collapsedFaces, const HashTable< label, label, Hash< label > > &edgeToEdge, const HashTable< label, label, Hash< label > > &edgeToFace, const label faceI, const label edgeI)
Calculates (cos of) angle across edgeI of faceI,.
Definition: triSurfaceTools.C:532
cut
Patchify triangles based on orientation w.r.t other (triangulated or triangulatable) surfaces.
Foam::PrimitivePatch::nPoints
label nPoints() const
Return number of points supporting patch faces.
Definition: PrimitivePatchTemplate.H:293
Foam::Vector::x
const Cmpt & x() const
Definition: VectorI.H:65
Foam::FatalError
error FatalError
Foam::triangle::b
const Point & b() const
Return second vertex.
Definition: triangleI.H:89
Foam::toPoint
PointFrompoint toPoint(const Foam::point &p)
Definition: pointConversion.H:80
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::triangle::POINT
@ POINT
Definition: triangle.H:97
Foam::geometryBase::name
const word & name() const
Return the name.
Definition: geometryBase.C:124
Foam::Vector::z
const Cmpt & z() const
Definition: VectorI.H:77
Foam::triSurfaceTools::NOEDGE
static const label NOEDGE
Definition: triSurfaceTools.H:367
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::PointIndexHit::setPoint
void setPoint(const Point &p)
Definition: PointIndexHit.H:163
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::triSurfaceTools::vertexUsesFace
static label vertexUsesFace(const triSurface &surf, const labelHashSet &faceUsed, const label vertI)
Definition: triSurfaceTools.C:427
Foam::triSurfaceTools::getVertexTriangles
static void getVertexTriangles(const triSurface &surf, const label edgeI, labelList &edgeTris)
Get all triangles using edge endpoint.
Definition: triSurfaceTools.C:1302
Foam::triSurfaceTools::visitFaces
static surfaceLocation visitFaces(const triSurface &s, const labelList &eFaces, const surfaceLocation &start, const label excludeEdgeI, const label excludePointI, const surfaceLocation &end, const plane &cutPlane)
Visits faces eFaces around start. Does not visit triangle.
Definition: triSurfaceTools.C:1166
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::OFstream
Output to file stream.
Definition: OFstream.H:81
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::triSurfaceTools::edgeSide
static sideType edgeSide(const triSurface &surf, const point &sample, const point &nearestPoint, const label edgeI)
If nearest point is on edgeI, determine on which side of surface.
Definition: triSurfaceTools.C:2141
Foam::triangle::nearestPointClassify
pointHit nearestPointClassify(const point &p, label &nearType, label &nearLabel) const
Find the nearest point to p on the triangle and classify it:
Definition: triangleI.H:521
Foam::PrimitivePatch::pointNormals
const Field< PointType > & pointNormals() const
Return point normals for patch.
Definition: PrimitivePatchTemplate.C:540
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::triSurfaceTools::trackToEdge
static surfaceLocation trackToEdge(const triSurface &, const surfaceLocation &start, const surfaceLocation &end, const plane &cutPlane)
Track on surface to get closer to point.
Definition: triSurfaceTools.C:2797
Foam::triSurfaceTools::classify
static surfaceLocation classify(const triSurface &, const label triI, const point &trianglePoint)
Test point on plane of triangle to see if on edge or point or inside.
Definition: triSurfaceTools.C:2753
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::triangle::EDGE
@ EDGE
Definition: triangle.H:98
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::PrimitivePatch::faceNormals
const Field< PointType > & faceNormals() const
Return face normals for patch.
Definition: PrimitivePatchTemplate.C:520
Foam::triSurfaceTools::minEdge
static label minEdge(const triSurface &surf, const labelList &edgeIndices)
Returns element in edgeIndices with minimum length.
Definition: triSurfaceTools.C:1978
f
labelList f(nPoints)
Foam::meshTools::otherFace
label otherFace(const primitiveMesh &, const label cellI, const label faceI, const label edgeI)
Return face on cell using edgeI but not faceI. Throws error.
Definition: meshTools.C:532
Foam::triSurfaceTools::getCollapsedFaces
static labelHashSet getCollapsedFaces(const triSurface &surf, label edgeI)
Faces to collapse because of edge collapse.
Definition: triSurfaceTools.C:374
Foam::triSurfaceTools::getVertexVertices
static labelList getVertexVertices(const triSurface &surf, const edge &e)
Get all vertices (local numbering) connected to vertices of edge.
Definition: triSurfaceTools.C:1345
Foam::Vector< scalar >
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
Foam::FixedList< scalar, 3 >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::triSurfaceTools::calcInterpolationWeights
static void calcInterpolationWeights(const triPointRef &, const point &, FixedList< scalar, 3 > &weights)
Calculate linear interpolation weights for point (guaranteed to be.
Definition: triSurfaceTools.C:2589
Foam::plane::refPoint
const point & refPoint() const
Return or return plane base point.
Definition: plane.C:255
Foam::surface
Definition: surface.H:55
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::labelledTri
Triangle with additional region number.
Definition: labelledTri.H:49
mergePoints.H
Merge points. See below.
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::surfaceLocation::triangle
label & triangle()
Definition: surfaceLocation.H:117
Foam::PrimitivePatch::localFaces
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatchTemplate.C:372
Foam::Vector::y
const Cmpt & y() const
Definition: VectorI.H:71
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::PrimitivePatch::faceFaces
const labelListList & faceFaces() const
Return face-face addressing.
Definition: PrimitivePatchTemplate.C:272
Foam::C
Graphite solid properties.
Definition: C.H:57
Foam::triSurfaceTools::sideType
sideType
On which side of surface.
Definition: triSurfaceTools.H:436
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
triSurfaceTools.H
Foam::mergePoints
label mergePoints(const UList< Type > &points, const scalar mergeTol, const bool verbose, labelList &pointMap, const Type &origin=Type::zero)
Sorts and merges points. All points closer than/equal mergeTol get merged.
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::triSurfaceTools::collapseMinCosAngle
static scalar collapseMinCosAngle(const triSurface &surf, const label v1, const point &pt, const labelHashSet &collapsedFaces, const HashTable< label, label, Hash< label > > &edgeToEdge, const HashTable< label, label, Hash< label > > &edgeToFace)
Calculate minimum (cos of) edge angle using addressing from.
Definition: triSurfaceTools.C:630
Foam::triSurfaceTools::cutEdge
static surfaceLocation cutEdge(const triSurface &s, const label triI, const label excludeEdgeI, const label excludePointI, const point &triPoint, const plane &cutPlane, const point &toPoint)
Finds the triangle edge/point cut by the plane between.
Definition: triSurfaceTools.C:843
Foam::triSurfaceTools::greenRefine
static void greenRefine(const triSurface &surf, const label faceI, const label edgeI, const label newPointI, DynamicList< labelledTri > &newFaces)
Definition: triSurfaceTools.C:91
Foam::triSurfaceTools::track
static void track(const triSurface &, const surfaceLocation &endInfo, const plane &cutPlane, surfaceLocation &hitInfo)
Track from edge to edge across surface. Uses trackToEdge.
Definition: triSurfaceTools.C:2875
normal
A normal distribution model.
Foam::triPointRef
triangle< point, const point & > triPointRef
Definition: triangle.H:75
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88
Foam::triSurfaceTools::otherEdges
static void otherEdges(const triSurface &surf, const label faceI, const label edgeI, label &e1, label &e2)
Get the two edges on faceI counterclockwise after edgeI.
Definition: triSurfaceTools.C:1458