1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2004-2010 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "triSurfaceTools.H"
28 #include "triSurface.H"
30 #include "mergePoints.H"
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 const Foam::label Foam::triSurfaceTools::ANYEDGE = -1;
39 const Foam::label Foam::triSurfaceTools::NOEDGE = -2;
40 const Foam::label Foam::triSurfaceTools::COLLAPSED = -3;
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
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",
52 // FaceI gets refined ('red'). Propagate edge refinement.
53 void Foam::triSurfaceTools::calcRefineStatus
55 const triSurface& surf,
57 List<refineType>& refine
60 if (refine[faceI] == RED)
62 // Already marked for refinement. Do nothing.
66 // Not marked or marked for 'green' refinement. Refine.
69 const labelList& myNeighbours = surf.faceFaces()[faceI];
71 forAll(myNeighbours, myNeighbourI)
73 label neighbourFaceI = myNeighbours[myNeighbourI];
75 if (refine[neighbourFaceI] == GREEN)
77 // Change to red refinement and propagate
78 calcRefineStatus(surf, neighbourFaceI, refine);
80 else if (refine[neighbourFaceI] == NONE)
82 refine[neighbourFaceI] = GREEN;
89 // Split faceI along edgeI at position newPointI
90 void Foam::triSurfaceTools::greenRefine
92 const triSurface& surf,
95 const label newPointI,
96 DynamicList<labelledTri>& newFaces
99 const labelledTri& f = surf.localFaces()[faceI];
100 const edge& e = surf.edges()[edgeI];
102 // Find index of edge in face.
104 label fp0 = findIndex(f, e[0]);
105 label fp1 = f.fcIndex(fp0);
106 label fp2 = f.fcIndex(fp1);
110 // Edge oriented like face
158 // Refine all triangles marked for refinement.
159 Foam::triSurface Foam::triSurfaceTools::doRefine
161 const triSurface& surf,
162 const List<refineType>& refineStatus
165 // Storage for new points. (start after old points)
166 DynamicList<point> newPoints(surf.nPoints());
167 forAll(surf.localPoints(), pointI)
169 newPoints.append(surf.localPoints()[pointI]);
171 label newVertI = surf.nPoints();
173 // Storage for new faces
174 DynamicList<labelledTri> newFaces(surf.size());
177 // Point index for midpoint on edge
178 labelList edgeMid(surf.nEdges(), -1);
180 forAll(refineStatus, faceI)
182 if (refineStatus[faceI] == RED)
184 // Create new vertices on all edges to be refined.
185 const labelList& fEdges = surf.faceEdges()[faceI];
189 label edgeI = fEdges[i];
191 if (edgeMid[edgeI] == -1)
193 const edge& e = surf.edges()[edgeI];
195 // Create new point on mid of edge
200 surf.localPoints()[e.start()]
201 + surf.localPoints()[e.end()]
204 edgeMid[edgeI] = newVertI++;
208 // Now we have new mid edge vertices for all edges on face
209 // so create triangles for RED rerfinement.
211 const edgeList& edges = surf.edges();
219 edges[fEdges[0]].commonVertex(edges[fEdges[1]]),
230 edges[fEdges[1]].commonVertex(edges[fEdges[2]]),
241 edges[fEdges[2]].commonVertex(edges[fEdges[0]]),
260 // Create triangles for GREEN refinement.
263 const label edgeI = fEdges[i];
265 label otherFaceI = otherFace(surf, faceI, edgeI);
267 if ((otherFaceI != -1) && (refineStatus[otherFaceI] == GREEN))
282 // Copy unmarked triangles since keep original vertex numbering.
283 forAll(refineStatus, faceI)
285 if (refineStatus[faceI] == NONE)
287 newFaces.append(surf.localFaces()[faceI]);
295 // Transfer DynamicLists to straight ones.
296 pointField allPoints;
297 allPoints.transfer(newPoints);
300 return triSurface(newFaces, surf.patches(), allPoints, true);
304 // Angle between two neighbouring triangles,
305 // angle between shared-edge end points and left and right face end points
306 Foam::scalar Foam::triSurfaceTools::faceCosAngle
314 const vector common(pEnd - pStart);
315 const vector base0(pLeft - pStart);
316 const vector base1(pRight - pStart);
318 vector n0(common ^ base0);
321 vector n1(base1 ^ common);
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
332 void Foam::triSurfaceTools::protectNeighbours
334 const triSurface& surf,
336 labelList& faceStatus
339 // const labelList& myFaces = surf.pointFaces()[vertI];
340 // forAll(myFaces, i)
342 // label faceI = myFaces[i];
344 // if ((faceStatus[faceI] == ANYEDGE) || (faceStatus[faceI] >= 0))
346 // faceStatus[faceI] = NOEDGE;
350 const labelList& myEdges = surf.pointEdges()[vertI];
353 const labelList& myFaces = surf.edgeFaces()[myEdges[i]];
355 forAll(myFaces, myFaceI)
357 label faceI = myFaces[myFaceI];
359 if ((faceStatus[faceI] == ANYEDGE) || (faceStatus[faceI] >= 0))
361 faceStatus[faceI] = NOEDGE;
369 // Edge collapse helper functions
372 // Get all faces that will get collapsed if edgeI collapses.
373 Foam::labelHashSet Foam::triSurfaceTools::getCollapsedFaces
375 const triSurface& surf,
379 const edge& e = surf.edges()[edgeI];
380 label v1 = e.start();
383 // Faces using edge will certainly get collapsed.
384 const labelList& myFaces = surf.edgeFaces()[edgeI];
386 labelHashSet facesToBeCollapsed(2*myFaces.size());
388 forAll(myFaces, myFaceI)
390 facesToBeCollapsed.insert(myFaces[myFaceI]);
393 // From faces using v1 check if they share an edge with faces
395 // - share edge: are part of 'splay' tree and will collapse if edge
397 const labelList& v1Faces = surf.pointFaces()[v1];
399 forAll(v1Faces, v1FaceI)
401 label face1I = v1Faces[v1FaceI];
403 label otherEdgeI = oppositeEdge(surf, face1I, v1);
405 // Step across edge to other face
406 label face2I = otherFace(surf, face1I, otherEdgeI);
410 // found face on other side of edge. Now check if uses v2.
411 if (oppositeVertex(surf, face2I, otherEdgeI) == v2)
413 // triangles face1I and face2I form splay tree and will
415 facesToBeCollapsed.insert(face1I);
416 facesToBeCollapsed.insert(face2I);
421 return facesToBeCollapsed;
425 // Return value of faceUsed for faces using vertI
426 Foam::label Foam::triSurfaceTools::vertexUsesFace
428 const triSurface& surf,
429 const labelHashSet& faceUsed,
433 const labelList& myFaces = surf.pointFaces()[vertI];
435 forAll(myFaces, myFaceI)
437 label face1I = myFaces[myFaceI];
439 if (faceUsed.found(face1I))
448 // Calculate new edge-face addressing resulting from edge collapse
449 void Foam::triSurfaceTools::getMergedEdges
451 const triSurface& surf,
453 const labelHashSet& collapsedFaces,
454 HashTable<label, label, Hash<label> >& edgeToEdge,
455 HashTable<label, label, Hash<label> >& edgeToFace
458 const edge& e = surf.edges()[edgeI];
459 label v1 = e.start();
462 const labelList& v1Faces = surf.pointFaces()[v1];
463 const labelList& v2Faces = surf.pointFaces()[v2];
465 // Mark all (non collapsed) faces using v2
466 labelHashSet v2FacesHash(v2Faces.size());
468 forAll(v2Faces, v2FaceI)
470 if (!collapsedFaces.found(v2Faces[v2FaceI]))
472 v2FacesHash.insert(v2Faces[v2FaceI]);
477 forAll(v1Faces, v1FaceI)
479 label face1I = v1Faces[v1FaceI];
481 if (collapsedFaces.found(face1I))
487 // Check if face1I uses a vertex connected to a face using v2
500 //Pout<< "Face:" << surf.localFaces()[face1I] << " other vertices:"
501 // << vert1I << ' ' << vert2I << endl;
503 // Check vert1, vert2 for usage by v2Face.
505 label commonVert = vert1I;
506 label face2I = vertexUsesFace(surf, v2FacesHash, commonVert);
510 face2I = vertexUsesFace(surf, v2FacesHash, commonVert);
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);
519 edgeToEdge.insert(edge1I, edge2I);
520 edgeToEdge.insert(edge2I, edge1I);
522 edgeToFace.insert(edge1I, face2I);
523 edgeToFace.insert(edge2I, face1I);
529 // Calculates (cos of) angle across edgeI of faceI,
530 // taking into account updated addressing (resulting from edge collapse)
531 Foam::scalar Foam::triSurfaceTools::edgeCosAngle
533 const triSurface& surf,
536 const labelHashSet& collapsedFaces,
537 const HashTable<label, label, Hash<label> >& edgeToEdge,
538 const HashTable<label, label, Hash<label> >& edgeToFace,
543 const pointField& localPoints = surf.localPoints();
545 label A = surf.edges()[edgeI].start();
546 label B = surf.edges()[edgeI].end();
547 label C = oppositeVertex(surf, faceI, edgeI);
553 if (edgeToEdge.found(edgeI))
555 // Use merged addressing
556 label edge2I = edgeToEdge[edgeI];
557 face2I = edgeToFace[edgeI];
559 D = oppositeVertex(surf, face2I, edge2I);
563 // Use normal edge-face addressing
564 face2I = otherFace(surf, faceI, edgeI);
566 if ((face2I != -1) && !collapsedFaces.found(face2I))
568 D = oppositeVertex(surf, face2I, edgeI);
578 cosAngle = faceCosAngle
588 cosAngle = faceCosAngle
598 cosAngle = faceCosAngle
608 cosAngle = faceCosAngle
618 FatalErrorIn("edgeCosAngle")
619 << "face " << faceI << " does not use vertex "
620 << v1 << " of collapsed edge" << abort(FatalError);
627 //- Calculate minimum (cos of) edge angle using addressing from collapsing
629 Foam::scalar Foam::triSurfaceTools::collapseMinCosAngle
631 const triSurface& surf,
634 const labelHashSet& collapsedFaces,
635 const HashTable<label, label, Hash<label> >& edgeToEdge,
636 const HashTable<label, label, Hash<label> >& edgeToFace
639 const labelList& v1Faces = surf.pointFaces()[v1];
643 forAll(v1Faces, v1FaceI)
645 label faceI = v1Faces[v1FaceI];
647 if (collapsedFaces.found(faceI))
652 const labelList& myEdges = surf.faceEdges()[faceI];
654 forAll(myEdges, myEdgeI)
656 label edgeI = myEdges[myEdgeI];
681 // Calculate max dihedral angle after collapsing edge to v1 (at pt).
682 // Note that all edges of all faces using v1 are affected.
683 bool Foam::triSurfaceTools::collapseCreatesFold
685 const triSurface& surf,
688 const labelHashSet& collapsedFaces,
689 const HashTable<label, label, Hash<label> >& edgeToEdge,
690 const HashTable<label, label, Hash<label> >& edgeToFace,
694 const labelList& v1Faces = surf.pointFaces()[v1];
696 forAll(v1Faces, v1FaceI)
698 label faceI = v1Faces[v1FaceI];
700 if (collapsedFaces.found(faceI))
705 const labelList& myEdges = surf.faceEdges()[faceI];
707 forAll(myEdges, myEdgeI)
709 label edgeI = myEdges[myEdgeI];
736 //// Return true if collapsing edgeI creates triangles on top of each other.
737 //// Is when the triangles neighbouring the collapsed one already share
739 //bool Foam::triSurfaceTools::collapseCreatesDuplicates
741 // const triSurface& surf,
742 // const label edgeI,
743 // const labelHashSet& collapsedFaces
746 //Pout<< "duplicate : edgeI:" << edgeI << surf.edges()[edgeI]
747 // << " collapsedFaces:" << collapsedFaces.toc() << endl;
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.
754 // HashTable<label, label, Hash<label> > neighbours;
756 // labelList collapsed = collapsedFaces.toc();
758 // forAll(collapsed, collapseI)
760 // const label faceI = collapsed[collapseI];
762 // const labelList& myEdges = surf.faceEdges()[faceI];
764 // Pout<< "collapsing faceI:" << faceI << " uses edges:" << myEdges
767 // forAll(myEdges, myEdgeI)
769 // const labelList& myFaces = surf.edgeFaces()[myEdges[myEdgeI]];
771 // Pout<< "Edge " << myEdges[myEdgeI] << " is used by faces "
772 // << myFaces << endl;
774 // if ((myEdges[myEdgeI] != edgeI) && (myFaces.size() == 2))
776 // // Get the other face
777 // label neighbourFaceI = myFaces[0];
779 // if (neighbourFaceI == faceI)
781 // neighbourFaceI = myFaces[1];
784 // // Is 'outside' face if not in collapsedFaces itself
785 // if (!collapsedFaces.found(neighbourFaceI))
787 // neighbours.insert(neighbourFaceI, myEdges[myEdgeI]);
793 // // Now neighbours contains first layer of triangles outside of
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
802 // // Check if any of the faces in neighbours share an edge. (n^2)
804 // labelList neighbourList = neighbours.toc();
806 //Pout<< "edgeI:" << edgeI << " neighbourList:" << neighbourList << endl;
809 // forAll(neighbourList, i)
811 // const labelList& faceIEdges = surf.faceEdges()[neighbourList[i]];
813 // for (label j = i+1; j < neighbourList.size(); i++)
815 // const labelList& faceJEdges = surf.faceEdges()[neighbourList[j]];
817 // // Check if faceI and faceJ share an edge
818 // forAll(faceIEdges, fI)
820 // forAll(faceJEdges, fJ)
822 // Pout<< " comparing " << faceIEdges[fI] << " to "
823 // << faceJEdges[fJ] << endl;
825 // if (faceIEdges[fI] == faceJEdges[fJ])
833 // Pout<< "Found no match. Returning false" << endl;
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
842 Foam::surfaceLocation Foam::triSurfaceTools::cutEdge
846 const label excludeEdgeI,
847 const label excludePointI,
849 const point& triPoint,
850 const plane& cutPlane,
854 const pointField& points = s.points();
855 const labelledTri& f = s[triI];
856 const labelList& fEdges = s.faceEdges()[triI];
858 // Get normal distance to planeN
859 FixedList<scalar, 3> d;
864 d[fp] = (points[f[fp]]-cutPlane.refPoint()) & cutPlane.normal();
868 // Normalise and truncate
873 if (mag(d[i]) < 1E-6)
879 // Return information
882 if (excludePointI != -1)
884 // Excluded point. Test only opposite edge.
886 label fp0 = findIndex(s.localFaces()[triI], excludePointI);
890 FatalErrorIn("cutEdge(..)") << "excludePointI:" << excludePointI
891 << " localF:" << s.localFaces()[triI] << abort(FatalError);
894 label fp1 = f.fcIndex(fp0);
895 label fp2 = f.fcIndex(fp1);
901 cut.setPoint(points[f[fp1]]);
902 cut.elementType() = triPointRef::POINT;
903 cut.setIndex(s.localFaces()[triI][fp1]);
905 else if (d[fp2] == 0.0)
908 cut.setPoint(points[f[fp2]]);
909 cut.elementType() = triPointRef::POINT;
910 cut.setIndex(s.localFaces()[triI][fp2]);
914 (d[fp1] < 0 && d[fp2] < 0)
915 || (d[fp1] > 0 && d[fp2] > 0)
918 // Both same sign. Not crossing edge at all.
919 // cut already set to miss().
926 (d[fp2]*points[f[fp1]] - d[fp1]*points[f[fp2]])
929 cut.elementType() = triPointRef::EDGE;
930 cut.setIndex(fEdges[fp1]);
935 // Find the two intersections
936 FixedList<surfaceLocation, 2> inters;
941 label fp1 = f.fcIndex(fp0);
947 FatalErrorIn("cutEdge(..)")
948 << "problem : triangle has three intersections." << nl
949 << "triangle:" << f.tri(points)
950 << " d:" << d << abort(FatalError);
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]);
960 (d[fp0] < 0 && d[fp1] > 0)
961 || (d[fp0] > 0 && d[fp1] < 0)
966 FatalErrorIn("cutEdge(..)")
967 << "problem : triangle has three intersections." << nl
968 << "triangle:" << f.tri(points)
969 << " d:" << d << abort(FatalError);
971 inters[interI].setHit();
972 inters[interI].setPoint
974 (d[fp0]*points[f[fp1]] - d[fp1]*points[f[fp0]])
977 inters[interI].elementType() = triPointRef::EDGE;
978 inters[interI].setIndex(fEdges[fp0]);
988 else if (interI == 1)
990 // Only one intersection. Should not happen!
993 else if (interI == 2)
995 // Handle excludeEdgeI
998 inters[0].elementType() == triPointRef::EDGE
999 && inters[0].index() == excludeEdgeI
1006 inters[1].elementType() == triPointRef::EDGE
1007 && inters[1].index() == excludeEdgeI
1014 // Two cuts. Find nearest.
1017 magSqr(inters[0].rawPoint() - toPoint)
1018 < magSqr(inters[1].rawPoint() - toPoint)
1034 // 'Snap' point on to endPoint.
1035 void Foam::triSurfaceTools::snapToEnd
1037 const triSurface& s,
1038 const surfaceLocation& end,
1039 surfaceLocation& current
1042 if (end.elementType() == triPointRef::NONE)
1044 if (current.elementType() == triPointRef::NONE)
1046 // endpoint on triangle; current on triangle
1047 if (current.index() == end.index())
1051 // Pout<< "snapToEnd : snapping:" << current << " onto:"
1058 // No need to handle current on edge/point since tracking handles this.
1060 else if (end.elementType() == triPointRef::EDGE)
1062 if (current.elementType() == triPointRef::NONE)
1064 // endpoint on edge; current on triangle
1065 const labelList& fEdges = s.faceEdges()[current.index()];
1067 if (findIndex(fEdges, end.index()) != -1)
1071 // Pout<< "snapToEnd : snapping:" << current << " onto:"
1078 else if (current.elementType() == triPointRef::EDGE)
1080 // endpoint on edge; current on edge
1081 if (current.index() == end.index())
1085 // Pout<< "snapToEnd : snapping:" << current << " onto:"
1094 // endpoint on edge; current on point
1095 const edge& e = s.edges()[end.index()];
1097 if (current.index() == e[0] || current.index() == e[1])
1101 // Pout<< "snapToEnd : snapping:" << current << " onto:"
1109 else // end.elementType() == POINT
1111 if (current.elementType() == triPointRef::NONE)
1113 // endpoint on point; current on triangle
1114 const triSurface::FaceType& f = s.localFaces()[current.index()];
1116 if (findIndex(f, end.index()) != -1)
1120 // Pout<< "snapToEnd : snapping:" << current << " onto:"
1127 else if (current.elementType() == triPointRef::EDGE)
1129 // endpoint on point; current on edge
1130 const edge& e = s.edges()[current.index()];
1132 if (end.index() == e[0] || end.index() == e[1])
1136 // Pout<< "snapToEnd : snapping:" << current << " onto:"
1145 // endpoint on point; current on point
1146 if (current.index() == end.index())
1150 // Pout<< "snapToEnd : snapping:" << current << " onto:"
1163 // - element type (triangle/edge/point) and index
1164 // - triangle to exclude
1165 Foam::surfaceLocation Foam::triSurfaceTools::visitFaces
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
1176 surfaceLocation nearest;
1178 scalar minDistSqr = Foam::sqr(GREAT);
1182 label triI = eFaces[i];
1184 // Make sure we don't revisit previous face
1185 if (triI != start.triangle())
1187 if (end.elementType() == triPointRef::NONE && end.index() == triI)
1189 // Endpoint is in this triangle. Jump there.
1192 nearest.triangle() = triI;
1197 // Which edge is cut.
1199 surfaceLocation cutInfo = cutEdge
1203 excludeEdgeI, // excludeEdgeI
1204 excludePointI, // excludePointI
1210 // If crossing an edge we expect next edge to be cut.
1211 if (excludeEdgeI != -1 && !cutInfo.hit())
1213 FatalErrorIn("triSurfaceTools::visitFaces(..)")
1214 << "Triangle:" << triI
1215 << " excludeEdge:" << excludeEdgeI
1216 << " point:" << start.rawPoint()
1217 << " plane:" << cutPlane
1218 << " . No intersection!" << abort(FatalError);
1223 scalar distSqr = magSqr(cutInfo.rawPoint()-end.rawPoint());
1225 if (distSqr < minDistSqr)
1227 minDistSqr = distSqr;
1229 nearest.triangle() = triI;
1237 if (nearest.triangle() == -1)
1239 // Did not move from edge. Give warning? Return something special?
1240 // For now responsability of caller to make sure that nothing has
1248 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1251 // Write pointField to file
1252 void Foam::triSurfaceTools::writeOBJ
1254 const fileName& fName,
1255 const pointField& pts
1258 OFstream outFile(fName);
1262 const point& pt = pts[pointI];
1264 outFile<< "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
1266 Pout<< "Written " << pts.size() << " vertices to file " << fName << endl;
1270 // Write vertex subset to OBJ format file
1271 void Foam::triSurfaceTools::writeOBJ
1273 const triSurface& surf,
1274 const fileName& fName,
1275 const boolList& markedVerts
1278 OFstream outFile(fName);
1281 forAll(markedVerts, vertI)
1283 if (markedVerts[vertI])
1285 const point& pt = surf.localPoints()[vertI];
1287 outFile<< "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
1292 Pout<< "Written " << nVerts << " vertices to file " << fName << endl;
1296 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1297 // Addressing helper functions:
1300 // Get all triangles using vertices of edge
1301 void Foam::triSurfaceTools::getVertexTriangles
1303 const triSurface& surf,
1308 const edge& e = surf.edges()[edgeI];
1309 const labelList& myFaces = surf.edgeFaces()[edgeI];
1311 label face1I = myFaces[0];
1313 if (myFaces.size() == 2)
1315 face2I = myFaces[1];
1318 const labelList& startFaces = surf.pointFaces()[e.start()];
1319 const labelList& endFaces = surf.pointFaces()[e.end()];
1321 // Number of triangles is sum of pointfaces - common faces
1322 // (= faces using edge)
1323 edgeTris.setSize(startFaces.size() + endFaces.size() - myFaces.size());
1326 forAll(startFaces, startFaceI)
1328 edgeTris[nTris++] = startFaces[startFaceI];
1331 forAll(endFaces, endFaceI)
1333 label faceI = endFaces[endFaceI];
1335 if ((faceI != face1I) && (faceI != face2I))
1337 edgeTris[nTris++] = faceI;
1343 // Get all vertices connected to vertices of edge
1344 Foam::labelList Foam::triSurfaceTools::getVertexVertices
1346 const triSurface& surf,
1350 const edgeList& edges = surf.edges();
1352 label v1 = e.start();
1355 // Get all vertices connected to v1 or v2 through an edge
1356 labelHashSet vertexNeighbours;
1358 const labelList& v1Edges = surf.pointEdges()[v1];
1360 forAll(v1Edges, v1EdgeI)
1362 const edge& e = edges[v1Edges[v1EdgeI]];
1363 vertexNeighbours.insert(e.otherVertex(v1));
1366 const labelList& v2Edges = surf.pointEdges()[v2];
1368 forAll(v2Edges, v2EdgeI)
1370 const edge& e = edges[v2Edges[v2EdgeI]];
1372 label vertI = e.otherVertex(v2);
1374 vertexNeighbours.insert(vertI);
1376 return vertexNeighbours.toc();
1380 //// Order vertices consistent with face
1381 //void Foam::triSurfaceTools::orderVertices
1383 // const labelledTri& f,
1390 // // Order v1, v2 in anticlockwise order.
1391 // bool reverse = false;
1400 // else if (f[1] == v1)
1428 // Get the other face using edgeI
1429 Foam::label Foam::triSurfaceTools::otherFace
1431 const triSurface& surf,
1436 const labelList& myFaces = surf.edgeFaces()[edgeI];
1438 if (myFaces.size() != 2)
1444 if (faceI == myFaces[0])
1456 // Get the two edges on faceI counterclockwise after edgeI
1457 void Foam::triSurfaceTools::otherEdges
1459 const triSurface& surf,
1466 const labelList& eFaces = surf.faceEdges()[faceI];
1468 label i0 = findIndex(eFaces, edgeI);
1475 "(const triSurface&, const label, const label,"
1477 ) << "Edge " << surf.edges()[edgeI] << " not in face "
1478 << surf.localFaces()[faceI] << abort(FatalError);
1481 label i1 = eFaces.fcIndex(i0);
1482 label i2 = eFaces.fcIndex(i1);
1489 // Get the two vertices on faceI counterclockwise vertI
1490 void Foam::triSurfaceTools::otherVertices
1492 const triSurface& surf,
1499 const labelledTri& f = surf.localFaces()[faceI];
1506 else if (vertI == f[1])
1511 else if (vertI == f[2])
1521 "(const triSurface&, const label, const label,"
1523 ) << "Vertex " << vertI << " not in face " << f << abort(FatalError);
1528 // Get edge opposite vertex
1529 Foam::label Foam::triSurfaceTools::oppositeEdge
1531 const triSurface& surf,
1536 const labelList& myEdges = surf.faceEdges()[faceI];
1538 forAll(myEdges, myEdgeI)
1540 label edgeI = myEdges[myEdgeI];
1542 const edge& e = surf.edges()[edgeI];
1544 if ((e.start() != vertI) && (e.end() != vertI))
1553 "(const triSurface&, const label, const label)"
1554 ) << "Cannot find vertex " << vertI << " in edges of face " << faceI
1555 << abort(FatalError);
1561 // Get vertex opposite edge
1562 Foam::label Foam::triSurfaceTools::oppositeVertex
1564 const triSurface& surf,
1569 const triSurface::FaceType& f = surf.localFaces()[faceI];
1570 const edge& e = surf.edges()[edgeI];
1574 label vertI = f[fp];
1576 if (vertI != e.start() && vertI != e.end())
1582 FatalErrorIn("triSurfaceTools::oppositeVertex")
1583 << "Cannot find vertex opposite edge " << edgeI << " vertices " << e
1584 << " in face " << faceI << " vertices " << f << abort(FatalError);
1590 // Returns edge label connecting v1, v2
1591 Foam::label Foam::triSurfaceTools::getEdge
1593 const triSurface& surf,
1598 const labelList& v1Edges = surf.pointEdges()[v1];
1600 forAll(v1Edges, v1EdgeI)
1602 label edgeI = v1Edges[v1EdgeI];
1603 const edge& e = surf.edges()[edgeI];
1605 if ((e.start() == v2) || (e.end() == v2))
1614 // Return index of triangle (or -1) using all three edges
1615 Foam::label Foam::triSurfaceTools::getTriangle
1617 const triSurface& surf,
1623 if ((e0I == e1I) || (e0I == e2I) || (e1I == e2I))
1628 "(const triSurface&, const label, const label,"
1630 ) << "Duplicate edge labels : e0:" << e0I << " e1:" << e1I
1632 << abort(FatalError);
1635 const labelList& eFaces = surf.edgeFaces()[e0I];
1637 forAll(eFaces, eFaceI)
1639 label faceI = eFaces[eFaceI];
1641 const labelList& myEdges = surf.faceEdges()[faceI];
1646 || (myEdges[1] == e1I)
1647 || (myEdges[2] == e1I)
1653 || (myEdges[1] == e2I)
1654 || (myEdges[2] == e2I)
1665 // Collapse indicated edges. Return new tri surface.
1666 Foam::triSurface Foam::triSurfaceTools::collapseEdges
1668 const triSurface& surf,
1669 const labelList& collapsableEdges
1672 pointField edgeMids(surf.nEdges());
1674 forAll(edgeMids, edgeI)
1676 const edge& e = surf.edges()[edgeI];
1681 surf.localPoints()[e.start()]
1682 + surf.localPoints()[e.end()]
1687 labelList faceStatus(surf.size(), ANYEDGE);
1689 //// Protect triangles which are on the border of different regions
1690 //forAll(edges, edgeI)
1692 // const labelList& neighbours = edgeFaces[edgeI];
1694 // if ((neighbours.size() != 2) && (neighbours.size() != 1))
1696 // FatalErrorIn("collapseEdges")
1697 // << abort(FatalError);
1700 // if (neighbours.size() == 2)
1702 // if (surf[neighbours[0]].region() != surf[neighbours[1]].region())
1704 // // Neighbours on different regions. For now dont allow
1706 // //Pout<< "protecting face " << neighbours[0]
1707 // // << ' ' << neighbours[1] << endl;
1708 // faceStatus[neighbours[0]] = NOEDGE;
1709 // faceStatus[neighbours[1]] = NOEDGE;
1714 return collapseEdges(surf, collapsableEdges, edgeMids, faceStatus);
1718 // Collapse indicated edges. Return new tri surface.
1719 Foam::triSurface Foam::triSurfaceTools::collapseEdges
1721 const triSurface& surf,
1722 const labelList& collapseEdgeLabels,
1723 const pointField& edgeMids,
1724 labelList& faceStatus
1727 const labelListList& edgeFaces = surf.edgeFaces();
1728 const pointField& localPoints = surf.localPoints();
1729 const edgeList& edges = surf.edges();
1731 // Storage for new points
1732 pointField newPoints(localPoints);
1734 // Map for old to new points
1735 labelList pointMap(localPoints.size());
1736 forAll(localPoints, pointI)
1738 pointMap[pointI] = pointI;
1742 // Do actual 'collapsing' of edges
1744 forAll(collapseEdgeLabels, collapseEdgeI)
1746 const label edgeI = collapseEdgeLabels[collapseEdgeI];
1748 if ((edgeI < 0) || (edgeI >= surf.nEdges()))
1750 FatalErrorIn("collapseEdges")
1751 << "Edge label outside valid range." << endl
1752 << "edge label:" << edgeI << endl
1753 << "total number of edges:" << surf.nEdges() << endl
1754 << abort(FatalError);
1757 const labelList& neighbours = edgeFaces[edgeI];
1759 if (neighbours.size() == 2)
1761 const label stat0 = faceStatus[neighbours[0]];
1762 const label stat1 = faceStatus[neighbours[1]];
1764 // Check faceStatus to make sure this one can be collapsed
1767 ((stat0 == ANYEDGE) || (stat0 == edgeI))
1768 && ((stat1 == ANYEDGE) || (stat1 == edgeI))
1771 const edge& e = edges[edgeI];
1773 // Set up mapping to 'collapse' points of edge
1776 (pointMap[e.start()] != e.start())
1777 || (pointMap[e.end()] != e.end())
1780 FatalErrorIn("collapseEdges")
1781 << "points already mapped. Double collapse." << endl
1782 << "edgeI:" << edgeI
1783 << " start:" << e.start()
1784 << " end:" << e.end()
1785 << " pointMap[start]:" << pointMap[e.start()]
1786 << " pointMap[end]:" << pointMap[e.end()]
1787 << abort(FatalError);
1790 const label minVert = min(e.start(), e.end());
1791 pointMap[e.start()] = minVert;
1792 pointMap[e.end()] = minVert;
1794 // Move shared vertex to mid of edge
1795 newPoints[minVert] = edgeMids[edgeI];
1797 // Protect neighbouring faces
1798 protectNeighbours(surf, e.start(), faceStatus);
1799 protectNeighbours(surf, e.end(), faceStatus);
1803 oppositeVertex(surf, neighbours[0], edgeI),
1809 oppositeVertex(surf, neighbours[1], edgeI),
1813 // Mark all collapsing faces
1814 labelList collapseFaces =
1821 forAll(collapseFaces, collapseI)
1823 faceStatus[collapseFaces[collapseI]] = COLLAPSED;
1830 // Storage for new triangles
1831 List<labelledTri> newTris(surf.size());
1834 const List<labelledTri>& localFaces = surf.localFaces();
1837 // Get only non-collapsed triangles and renumber vertex labels.
1838 forAll(localFaces, faceI)
1840 const labelledTri& f = localFaces[faceI];
1842 const label a = pointMap[f[0]];
1843 const label b = pointMap[f[1]];
1844 const label c = pointMap[f[2]];
1848 (a != b) && (a != c) && (b != c)
1849 && (faceStatus[faceI] != COLLAPSED)
1852 // uncollapsed triangle
1853 newTris[newTriI++] = labelledTri(a, b, c, f.region());
1857 //Pout<< "Collapsed triangle " << faceI
1858 // << " vertices:" << f << endl;
1861 newTris.setSize(newTriI);
1867 triSurface tempSurf(newTris, surf.patches(), newPoints);
1872 tempSurf.localFaces(),
1874 tempSurf.localPoints()
1879 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1881 Foam::triSurface Foam::triSurfaceTools::redGreenRefine
1883 const triSurface& surf,
1884 const labelList& refineFaces
1887 List<refineType> refineStatus(surf.size(), NONE);
1889 // Mark & propagate refinement
1890 forAll(refineFaces, refineFaceI)
1892 calcRefineStatus(surf, refineFaces[refineFaceI], refineStatus);
1895 // Do actual refinement
1896 return doRefine(surf, refineStatus);
1900 Foam::triSurface Foam::triSurfaceTools::greenRefine
1902 const triSurface& surf,
1903 const labelList& refineEdges
1906 // Storage for marking faces
1907 List<refineType> refineStatus(surf.size(), NONE);
1909 // Storage for new faces
1910 DynamicList<labelledTri> newFaces(0);
1912 pointField newPoints(surf.localPoints());
1913 newPoints.setSize(surf.nPoints() + surf.nEdges());
1914 label newPointI = surf.nPoints();
1918 forAll(refineEdges, refineEdgeI)
1920 label edgeI = refineEdges[refineEdgeI];
1922 const labelList& myFaces = surf.edgeFaces()[edgeI];
1924 bool neighbourIsRefined= false;
1926 forAll(myFaces, myFaceI)
1928 if (refineStatus[myFaces[myFaceI]] != NONE)
1930 neighbourIsRefined = true;
1934 // Only refine if none of the faces is refined
1935 if (!neighbourIsRefined)
1938 const edge& e = surf.edges()[edgeI];
1943 surf.localPoints()[e.start()]
1944 + surf.localPoints()[e.end()]
1947 newPoints[newPointI] = mid;
1949 // Refine faces using edge
1950 forAll(myFaces, myFaceI)
1952 // Add faces to newFaces
1963 refineStatus[myFaces[myFaceI]] = GREEN;
1970 // Add unrefined faces
1971 forAll(surf.localFaces(), faceI)
1973 if (refineStatus[faceI] == NONE)
1975 newFaces.append(surf.localFaces()[faceI]);
1980 newPoints.setSize(newPointI);
1982 return triSurface(newFaces, surf.patches(), newPoints, true);
1987 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1988 // Geometric helper functions:
1991 // Returns element in edgeIndices with minimum length
1992 Foam::label Foam::triSurfaceTools::minEdge
1994 const triSurface& surf,
1995 const labelList& edgeIndices
1998 scalar minLength = GREAT;
1999 label minIndex = -1;
2000 forAll(edgeIndices, i)
2002 const edge& e = surf.edges()[edgeIndices[i]];
2007 surf.localPoints()[e.end()]
2008 - surf.localPoints()[e.start()]
2011 if (length < minLength)
2017 return edgeIndices[minIndex];
2021 // Returns element in edgeIndices with maximum length
2022 Foam::label Foam::triSurfaceTools::maxEdge
2024 const triSurface& surf,
2025 const labelList& edgeIndices
2028 scalar maxLength = -GREAT;
2029 label maxIndex = -1;
2030 forAll(edgeIndices, i)
2032 const edge& e = surf.edges()[edgeIndices[i]];
2037 surf.localPoints()[e.end()]
2038 - surf.localPoints()[e.start()]
2041 if (length > maxLength)
2047 return edgeIndices[maxIndex];
2051 // Merge points and reconstruct surface
2052 Foam::triSurface Foam::triSurfaceTools::mergePoints
2054 const triSurface& surf,
2055 const scalar mergeTol
2058 pointField newPoints(surf.nPoints());
2060 labelList pointMap(surf.nPoints());
2062 bool hasMerged = Foam::mergePoints
2073 // Pack the triangles
2075 // Storage for new triangles
2076 List<labelledTri> newTriangles(surf.size());
2077 label newTriangleI = 0;
2081 const labelledTri& f = surf.localFaces()[faceI];
2083 label newA = pointMap[f[0]];
2084 label newB = pointMap[f[1]];
2085 label newC = pointMap[f[2]];
2087 if ((newA != newB) && (newA != newC) && (newB != newC))
2089 newTriangles[newTriangleI++] =
2090 labelledTri(newA, newB, newC, f.region());
2093 newTriangles.setSize(newTriangleI);
2109 // Calculate normal on triangle
2110 Foam::vector Foam::triSurfaceTools::surfaceNormal
2112 const triSurface& surf,
2113 const label nearestFaceI,
2114 const point& nearestPt
2117 const triSurface::FaceType& f = surf[nearestFaceI];
2118 const pointField& points = surf.points();
2120 label nearType, nearLabel;
2122 f.nearestPointClassify(nearestPt, points, nearType, nearLabel);
2124 if (nearType == triPointRef::NONE)
2127 return surf.faceNormals()[nearestFaceI];
2129 else if (nearType == triPointRef::EDGE)
2131 // Nearest to edge. Assume order of faceEdges same as face vertices.
2132 label edgeI = surf.faceEdges()[nearestFaceI][nearLabel];
2134 // Calculate edge normal by averaging face normals
2135 const labelList& eFaces = surf.edgeFaces()[edgeI];
2137 vector edgeNormal(vector::zero);
2141 edgeNormal += surf.faceNormals()[eFaces[i]];
2143 return edgeNormal/(mag(edgeNormal) + VSMALL);
2148 const triSurface::FaceType& localF = surf.localFaces()[nearestFaceI];
2149 return surf.pointNormals()[localF[nearLabel]];
2154 Foam::triSurfaceTools::sideType Foam::triSurfaceTools::edgeSide
2156 const triSurface& surf,
2157 const point& sample,
2158 const point& nearestPoint,
2162 const labelList& eFaces = surf.edgeFaces()[edgeI];
2164 if (eFaces.size() != 2)
2166 // Surface not closed.
2171 const vectorField& faceNormals = surf.faceNormals();
2173 // Compare to bisector. This is actually correct since edge is
2174 // nearest so there is a knife-edge.
2176 vector n = 0.5*(faceNormals[eFaces[0]] + faceNormals[eFaces[1]]);
2178 if (((sample - nearestPoint) & n) > 0)
2190 // Calculate normal on triangle
2191 Foam::triSurfaceTools::sideType Foam::triSurfaceTools::surfaceSide
2193 const triSurface& surf,
2194 const point& sample,
2195 const label nearestFaceI
2198 const triSurface::FaceType& f = surf[nearestFaceI];
2199 const pointField& points = surf.points();
2201 // Find where point is on face
2202 label nearType, nearLabel;
2204 pointHit pHit = f.nearestPointClassify(sample, points, nearType, nearLabel);
2206 const point& nearestPoint(pHit.rawPoint());
2208 if (nearType == triPointRef::NONE)
2210 vector sampleNearestVec = (sample - nearestPoint);
2212 // Nearest to face interior. Use faceNormal to determine side
2213 scalar c = sampleNearestVec & surf.faceNormals()[nearestFaceI];
2215 // // If the sample is essentially on the face, do not check for
2216 // // it being perpendicular.
2218 // scalar magSampleNearestVec = mag(sampleNearestVec);
2220 // if (magSampleNearestVec > SMALL)
2222 // c /= magSampleNearestVec*mag(surf.faceNormals()[nearestFaceI]);
2224 // if (mag(c) < 0.99)
2226 // FatalErrorIn("triSurfaceTools::surfaceSide")
2227 // << "nearestPoint identified as being on triangle face "
2228 // << "but vector from nearestPoint to sample is not "
2229 // << "perpendicular to the normal." << nl
2230 // << "sample: " << sample << nl
2231 // << "nearestPoint: " << nearestPoint << nl
2232 // << "sample - nearestPoint: "
2233 // << sample - nearestPoint << nl
2234 // << "normal: " << surf.faceNormals()[nearestFaceI] << nl
2235 // << "mag(sample - nearestPoint): "
2236 // << mag(sample - nearestPoint) << nl
2237 // << "normalised dot product: " << c << nl
2238 // << "triangle vertices: " << nl
2239 // << " " << points[f[0]] << nl
2240 // << " " << points[f[1]] << nl
2241 // << " " << points[f[2]] << nl
2242 // << abort(FatalError);
2255 else if (nearType == triPointRef::EDGE)
2257 // Nearest to edge nearLabel. Note that this can only be a knife-edge
2258 // situation since otherwise the nearest point could never be the edge.
2260 // Get the edge. Assume order of faceEdges same as face vertices.
2261 label edgeI = surf.faceEdges()[nearestFaceI][nearLabel];
2265 // // Check order of faceEdges same as face vertices.
2266 // const edge& e = surf.edges()[edgeI];
2267 // const labelList& meshPoints = surf.meshPoints();
2268 // const edge meshEdge(meshPoints[e[0]], meshPoints[e[1]]);
2273 // != edge(f[nearLabel], f[f.fcIndex(nearLabel)])
2276 // FatalErrorIn("triSurfaceTools::surfaceSide")
2277 // << "Edge:" << edgeI << " local vertices:" << e
2278 // << " mesh vertices:" << meshEdge
2279 // << " not at position " << nearLabel
2280 // << " in face " << f
2281 // << abort(FatalError);
2285 return edgeSide(surf, sample, nearestPoint, edgeI);
2289 // Nearest to point. Could use pointNormal here but is not correct.
2290 // Instead determine which edge using point is nearest and use test
2291 // above (nearType == triPointRef::EDGE).
2294 const triSurface::FaceType& localF = surf.localFaces()[nearestFaceI];
2295 label nearPointI = localF[nearLabel];
2297 const edgeList& edges = surf.edges();
2298 const pointField& localPoints = surf.localPoints();
2299 const point& base = localPoints[nearPointI];
2301 const labelList& pEdges = surf.pointEdges()[nearPointI];
2303 scalar minDistSqr = Foam::sqr(GREAT);
2304 label minEdgeI = -1;
2308 label edgeI = pEdges[i];
2310 const edge& e = edges[edgeI];
2312 label otherPointI = e.otherVertex(nearPointI);
2315 vector eVec(localPoints[otherPointI] - base);
2316 scalar magEVec = mag(eVec);
2318 if (magEVec > VSMALL)
2322 // Get point along vector and determine closest.
2323 const point perturbPoint = base + eVec;
2325 scalar distSqr = Foam::magSqr(sample - perturbPoint);
2327 if (distSqr < minDistSqr)
2329 minDistSqr = distSqr;
2337 FatalErrorIn("treeDataTriSurface::getSide")
2338 << "Problem: did not find edge closer than " << minDistSqr
2339 << abort(FatalError);
2342 return edgeSide(surf, sample, nearestPoint, minEdgeI);
2347 // triangulation of boundaryMesh
2348 Foam::triSurface Foam::triSurfaceTools::triangulate
2350 const polyBoundaryMesh& bMesh,
2351 const labelHashSet& includePatches,
2355 const polyMesh& mesh = bMesh.mesh();
2357 // Storage for surfaceMesh. Size estimate.
2358 DynamicList<labelledTri> triangles
2360 mesh.nFaces() - mesh.nInternalFaces()
2363 label newPatchI = 0;
2365 forAllConstIter(labelHashSet, includePatches, iter)
2367 const label patchI = iter.key();
2368 const polyPatch& patch = bMesh[patchI];
2369 const pointField& points = patch.points();
2371 label nTriTotal = 0;
2373 forAll(patch, patchFaceI)
2375 const face& f = patch[patchFaceI];
2377 faceList triFaces(f.nTriangles(points));
2381 f.triangles(points, nTri, triFaces);
2383 forAll(triFaces, triFaceI)
2385 const face& f = triFaces[triFaceI];
2387 triangles.append(labelledTri(f[0], f[1], f[2], newPatchI));
2395 Pout<< patch.name() << " : generated " << nTriTotal
2396 << " triangles from " << patch.size() << " faces with"
2397 << " new patchid " << newPatchI << endl;
2404 // Create globally numbered tri surface
2405 triSurface rawSurface(triangles, mesh.points());
2407 // Create locally numbered tri surface
2410 rawSurface.localFaces(),
2411 rawSurface.localPoints()
2414 // Add patch names to surface
2415 surface.patches().setSize(newPatchI);
2419 forAllConstIter(labelHashSet, includePatches, iter)
2421 const label patchI = iter.key();
2422 const polyPatch& patch = bMesh[patchI];
2424 surface.patches()[newPatchI].name() = patch.name();
2425 surface.patches()[newPatchI].geometricType() = patch.type();
2434 // triangulation of boundaryMesh
2435 Foam::triSurface Foam::triSurfaceTools::triangulateFaceCentre
2437 const polyBoundaryMesh& bMesh,
2438 const labelHashSet& includePatches,
2442 const polyMesh& mesh = bMesh.mesh();
2444 // Storage for new points = meshpoints + face centres.
2445 const pointField& points = mesh.points();
2446 const pointField& faceCentres = mesh.faceCentres();
2448 pointField newPoints(points.size() + faceCentres.size());
2450 label newPointI = 0;
2452 forAll(points, pointI)
2454 newPoints[newPointI++] = points[pointI];
2456 forAll(faceCentres, faceI)
2458 newPoints[newPointI++] = faceCentres[faceI];
2462 // Count number of faces.
2463 DynamicList<labelledTri> triangles
2465 mesh.nFaces() - mesh.nInternalFaces()
2468 label newPatchI = 0;
2470 forAllConstIter(labelHashSet, includePatches, iter)
2472 const label patchI = iter.key();
2473 const polyPatch& patch = bMesh[patchI];
2475 label nTriTotal = 0;
2477 forAll(patch, patchFaceI)
2479 // Face in global coords.
2480 const face& f = patch[patchFaceI];
2482 // Index in newPointI of face centre.
2483 label fc = points.size() + patchFaceI + patch.start();
2487 label fp1 = f.fcIndex(fp);
2489 triangles.append(labelledTri(f[fp], f[fp1], fc, newPatchI));
2497 Pout<< patch.name() << " : generated " << nTriTotal
2498 << " triangles from " << patch.size() << " faces with"
2499 << " new patchid " << newPatchI << endl;
2507 // Create globally numbered tri surface
2508 triSurface rawSurface(triangles, newPoints);
2510 // Create locally numbered tri surface
2513 rawSurface.localFaces(),
2514 rawSurface.localPoints()
2517 // Add patch names to surface
2518 surface.patches().setSize(newPatchI);
2522 forAllConstIter(labelHashSet, includePatches, iter)
2524 const label patchI = iter.key();
2525 const polyPatch& patch = bMesh[patchI];
2527 surface.patches()[newPatchI].name() = patch.name();
2528 surface.patches()[newPatchI].geometricType() = patch.type();
2537 Foam::triSurface Foam::triSurfaceTools::delaunay2D(const List<vector2D>& pts)
2539 // Vertices in geompack notation. Note that could probably just use
2540 // pts.begin() if double precision.
2541 List<doubleScalar> geompackVertices(2*pts.size());
2545 geompackVertices[doubleI++] = pts[i][0];
2546 geompackVertices[doubleI++] = pts[i][1];
2549 // Storage for triangles
2551 List<int> triangle_node(m2*3*pts.size());
2552 List<int> triangle_neighbor(m2*3*pts.size());
2559 geompackVertices.begin(),
2561 triangle_node.begin(),
2562 triangle_neighbor.begin()
2567 FatalErrorIn("triSurfaceTools::delaunay2D(const List<vector2D>&)")
2568 << "Failed dtris2 with vertices:" << pts.size()
2569 << abort(FatalError);
2573 triangle_node.setSize(3*nTris);
2574 triangle_neighbor.setSize(3*nTris);
2576 // Convert to triSurface.
2577 List<labelledTri> faces(nTris);
2581 faces[i] = labelledTri
2583 triangle_node[3*i]-1,
2584 triangle_node[3*i+1]-1,
2585 triangle_node[3*i+2]-1,
2590 pointField points(pts.size());
2593 points[i][0] = pts[i][0];
2594 points[i][1] = pts[i][1];
2598 return triSurface(faces, points);
2602 //// Use CGAL to do Delaunay
2603 //#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
2604 //#include <CGAL/Delaunay_triangulation_2.h>
2605 //#include <CGAL/Triangulation_vertex_base_with_info_2.h>
2606 //#include <cassert>
2608 //struct K : CGAL::Exact_predicates_inexact_constructions_kernel {};
2610 //typedef CGAL::Triangulation_vertex_base_with_info_2<Foam::label, K> Vb;
2611 //typedef CGAL::Triangulation_data_structure_2<Vb> Tds;
2612 //typedef CGAL::Delaunay_triangulation_2<K,Tds> Triangulation;
2614 //typedef Triangulation::Vertex_handle Vertex_handle;
2615 //typedef Triangulation::Vertex_iterator Vertex_iterator;
2616 //typedef Triangulation::Face_handle Face_handle;
2617 //typedef Triangulation::Finite_faces_iterator Finite_faces_iterator;
2618 //typedef Triangulation::Point Point;
2619 //Foam::triSurface Foam::triSurfaceTools::delaunay2D(const List<vector2D>& pts)
2623 // // Insert 2D vertices; building triangulation
2626 // const point& pt = pts[i];
2628 // T.insert(Point(pt[0], pt[1]));
2632 // // Number vertices
2633 // // ~~~~~~~~~~~~~~~
2638 // Vertex_iterator it = T.vertices_begin();
2639 // it != T.vertices_end();
2643 // it->info() = vertI++;
2649 // List<labelledTri> faces(T.number_of_faces());
2655 // Finite_faces_iterator fc = T.finite_faces_begin();
2656 // fc != T.finite_faces_end();
2660 // faces[faceI++] = labelledTri
2662 // fc->vertex(0)->info(),
2663 // f[1] = fc->vertex(1)->info(),
2664 // f[2] = fc->vertex(2)->info()
2668 // pointField points(pts.size());
2671 // points[i][0] = pts[i][0];
2672 // points[i][1] = pts[i][1];
2673 // points[i][2] = 0.0;
2676 // return triSurface(faces, points);
2680 void Foam::triSurfaceTools::calcInterpolationWeights
2682 const triPointRef& tri,
2684 FixedList<scalar, 3>& weights
2687 // calculate triangle edge vectors and triangle face normal
2688 // the 'i':th edge is opposite node i
2689 FixedList<vector, 3> edge;
2690 edge[0] = tri.c()-tri.b();
2691 edge[1] = tri.a()-tri.c();
2692 edge[2] = tri.b()-tri.a();
2694 vector triangleFaceNormal = edge[1] ^ edge[2];
2696 // calculate edge normal (pointing inwards)
2697 FixedList<vector, 3> normal;
2698 for (label i=0; i<3; i++)
2700 normal[i] = triangleFaceNormal ^ edge[i];
2701 normal[i] /= mag(normal[i]) + VSMALL;
2704 weights[0] = ((p-tri.b()) & normal[0]) / max(VSMALL, normal[0] & edge[1]);
2705 weights[1] = ((p-tri.c()) & normal[1]) / max(VSMALL, normal[1] & edge[2]);
2706 weights[2] = ((p-tri.a()) & normal[2]) / max(VSMALL, normal[2] & edge[0]);
2710 // Calculate weighting factors from samplePts to triangle it is in.
2711 // Uses linear search.
2712 void Foam::triSurfaceTools::calcInterpolationWeights
2714 const triSurface& s,
2715 const pointField& samplePts,
2716 List<FixedList<label, 3> >& allVerts,
2717 List<FixedList<scalar, 3> >& allWeights
2720 allVerts.setSize(samplePts.size());
2721 allWeights.setSize(samplePts.size());
2723 const pointField& points = s.points();
2725 forAll(samplePts, i)
2727 const point& samplePt = samplePts[i];
2730 FixedList<label, 3>& verts = allVerts[i];
2731 FixedList<scalar, 3>& weights = allWeights[i];
2733 scalar minDistance = GREAT;
2737 const labelledTri& f = s[faceI];
2739 triPointRef tri(f.tri(points));
2741 label nearType, nearLabel;
2743 pointHit nearest = tri.nearestPointClassify
2752 // samplePt inside triangle
2757 calcInterpolationWeights(tri, nearest.rawPoint(), weights);
2759 //Pout<< "calcScalingFactors : samplePt:" << samplePt
2760 // << " inside triangle:" << faceI
2761 // << " verts:" << verts
2762 // << " weights:" << weights
2767 else if (nearest.distance() < minDistance)
2769 minDistance = nearest.distance();
2771 // Outside triangle. Store nearest.
2773 if (nearType == triPointRef::POINT)
2775 verts[0] = f[nearLabel];
2778 weights[1] = -GREAT;
2780 weights[2] = -GREAT;
2782 //Pout<< "calcScalingFactors : samplePt:" << samplePt
2783 // << " distance:" << nearest.distance()
2784 // << " from point:" << points[f[nearLabel]]
2787 else if (nearType == triPointRef::EDGE)
2789 verts[0] = f[nearLabel];
2790 verts[1] = f[f.fcIndex(nearLabel)];
2793 const point& p0 = points[verts[0]];
2794 const point& p1 = points[verts[1]];
2802 mag(nearest.rawPoint() - p0)/mag(p1 - p0)
2809 weights[2] = -GREAT;
2811 //Pout<< "calcScalingFactors : samplePt:" << samplePt
2812 // << " distance:" << nearest.distance()
2813 // << " from edge:" << p0 << p1 << " s:" << s
2818 // triangle. Can only happen because of truncation errors.
2823 calcInterpolationWeights(tri, nearest.rawPoint(), weights);
2833 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2837 // Test point on surface to see if is on face,edge or point.
2838 Foam::surfaceLocation Foam::triSurfaceTools::classify
2840 const triSurface& s,
2842 const point& trianglePoint
2845 surfaceLocation nearest;
2847 // Nearest point could be on point or edge. Retest.
2848 label index, elemType;
2850 triPointRef(s[triI].tri(s.points())).classify
2857 nearest.setPoint(trianglePoint);
2859 if (elemType == triPointRef::NONE)
2862 nearest.setIndex(triI);
2863 nearest.elementType() = triPointRef::NONE;
2865 else if (elemType == triPointRef::EDGE)
2868 nearest.setIndex(s.faceEdges()[triI][index]);
2869 nearest.elementType() = triPointRef::EDGE;
2871 else //if (elemType == triPointRef::POINT)
2874 nearest.setIndex(s.localFaces()[triI][index]);
2875 nearest.elementType() = triPointRef::POINT;
2882 Foam::surfaceLocation Foam::triSurfaceTools::trackToEdge
2884 const triSurface& s,
2885 const surfaceLocation& start,
2886 const surfaceLocation& end,
2887 const plane& cutPlane
2890 // Start off from starting point
2891 surfaceLocation nearest = start;
2894 // See if in same triangle as endpoint. If so snap.
2895 snapToEnd(s, end, nearest);
2899 // Not yet at end point
2901 if (start.elementType() == triPointRef::NONE)
2903 // Start point is inside triangle. Trivial cases already handled
2906 // end point is on edge or point so cross currrent triangle to
2907 // see which edge is cut.
2912 start.index(), // triangle
2919 nearest.elementType() = triPointRef::EDGE;
2920 nearest.triangle() = start.index();
2923 else if (start.elementType() == triPointRef::EDGE)
2925 // Pick connected triangle that is most in direction.
2926 const labelList& eFaces = s.edgeFaces()[start.index()];
2928 nearest = visitFaces
2933 start.index(), // excludeEdgeI
2934 -1, // excludePointI
2939 else // start.elementType() == triPointRef::POINT
2941 const labelList& pFaces = s.pointFaces()[start.index()];
2943 nearest = visitFaces
2949 start.index(), // excludePointI
2954 snapToEnd(s, end, nearest);
2960 void Foam::triSurfaceTools::track
2962 const triSurface& s,
2963 const surfaceLocation& endInfo,
2964 const plane& cutPlane,
2965 surfaceLocation& hitInfo
2968 //OFstream str("track.obj");
2970 //meshTools::writeOBJ(str, hitInfo.rawPoint());
2973 // Track across surface.
2976 //Pout<< "Tracking from:" << nl
2977 // << " " << hitInfo.info()
2980 hitInfo = trackToEdge
2988 //meshTools::writeOBJ(str, hitInfo.rawPoint());
2990 //str<< "l " << vertI-1 << ' ' << vertI << nl;
2992 //Pout<< "Tracked to:" << nl
2993 // << " " << hitInfo.info() << endl;
2995 if (hitInfo.hit() || hitInfo.triangle() == -1)
3003 // ************************************************************************* //