Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / meshTools / triSurface / triSurfaceTools / triSurfaceTools.C
blobd9b29a3fff9718cc38fa288fc93a3bd585016c55
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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"
29 #include "OFstream.H"
30 #include "mergePoints.H"
31 #include "polyMesh.H"
32 #include "plane.H"
33 #include "geompack.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",
49     Wiley-Teubner, 1996)
52 // FaceI gets refined ('red'). Propagate edge refinement.
53 void Foam::triSurfaceTools::calcRefineStatus
55     const triSurface& surf,
56     const label faceI,
57     List<refineType>& refine
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;
69         const labelList& myNeighbours = surf.faceFaces()[faceI];
71         forAll(myNeighbours, myNeighbourI)
72         {
73             label neighbourFaceI = myNeighbours[myNeighbourI];
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     }
89 // Split faceI along edgeI at position newPointI
90 void Foam::triSurfaceTools::greenRefine
92     const triSurface& surf,
93     const label faceI,
94     const label edgeI,
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);
108     if (f[fp1] == e[1])
109     {
110         // Edge oriented like face
111         newFaces.append
112         (
113             labelledTri
114             (
115                 f[fp0],
116                 newPointI,
117                 f[fp2],
118                 f.region()
119             )
120         );
121         newFaces.append
122         (
123             labelledTri
124             (
125                 newPointI,
126                 f[fp1],
127                 f[fp2],
128                 f.region()
129             )
130         );
131     }
132     else
133     {
134         newFaces.append
135         (
136             labelledTri
137             (
138                 f[fp2],
139                 newPointI,
140                 f[fp1],
141                 f.region()
142             )
143         );
144         newFaces.append
145         (
146             labelledTri
147             (
148                 newPointI,
149                 f[fp0],
150                 f[fp1],
151                 f.region()
152             )
153         );
154     }
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)
168     {
169         newPoints.append(surf.localPoints()[pointI]);
170     }
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)
181     {
182         if (refineStatus[faceI] == RED)
183         {
184             // Create new vertices on all edges to be refined.
185             const labelList& fEdges = surf.faceEdges()[faceI];
187             forAll(fEdges, i)
188             {
189                 label edgeI = fEdges[i];
191                 if (edgeMid[edgeI] == -1)
192                 {
193                     const edge& e = surf.edges()[edgeI];
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             }
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();
213             // Corner triangles
214             newFaces.append
215             (
216                 labelledTri
217                 (
218                     edgeMid[fEdges[0]],
219                     edges[fEdges[0]].commonVertex(edges[fEdges[1]]),
220                     edgeMid[fEdges[1]],
221                     surf[faceI].region()
222                 )
223             );
225             newFaces.append
226             (
227                 labelledTri
228                 (
229                     edgeMid[fEdges[1]],
230                     edges[fEdges[1]].commonVertex(edges[fEdges[2]]),
231                     edgeMid[fEdges[2]],
232                     surf[faceI].region()
233                 )
234             );
236             newFaces.append
237             (
238                 labelledTri
239                 (
240                     edgeMid[fEdges[2]],
241                     edges[fEdges[2]].commonVertex(edges[fEdges[0]]),
242                     edgeMid[fEdges[0]],
243                     surf[faceI].region()
244                 )
245             );
247             // Inner triangle
248             newFaces.append
249             (
250                 labelledTri
251                 (
252                     edgeMid[fEdges[0]],
253                     edgeMid[fEdges[1]],
254                     edgeMid[fEdges[2]],
255                     surf[faceI].region()
256                 )
257             );
260             // Create triangles for GREEN refinement.
261             forAll(fEdges, i)
262             {
263                 const label edgeI = fEdges[i];
265                 label otherFaceI = otherFace(surf, faceI, edgeI);
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     }
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     }
291     newFaces.shrink();
292     newPoints.shrink();
295     // Transfer DynamicLists to straight ones.
296     pointField allPoints;
297     allPoints.transfer(newPoints);
298     newPoints.clear();
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
308     const point& pStart,
309     const point& pEnd,
310     const point& pLeft,
311     const point& pRight
314     const vector common(pEnd - pStart);
315     const vector base0(pLeft - pStart);
316     const vector base1(pRight - pStart);
318     vector n0(common ^ base0);
319     n0 /= Foam::mag(n0);
321     vector n1(base1 ^ common);
322     n1 /= Foam::mag(n1);
324     return n0 & n1;
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,
335     const label vertI,
336     labelList& faceStatus
339 //    const labelList& myFaces = surf.pointFaces()[vertI];
340 //    forAll(myFaces, i)
341 //    {
342 //        label faceI = myFaces[i];
344 //        if ((faceStatus[faceI] == ANYEDGE) || (faceStatus[faceI] >= 0))
345 //        {
346 //            faceStatus[faceI] = NOEDGE;
347 //        }
348 //    }
350     const labelList& myEdges = surf.pointEdges()[vertI];
351     forAll(myEdges, i)
352     {
353         const labelList& myFaces = surf.edgeFaces()[myEdges[i]];
355         forAll(myFaces, myFaceI)
356         {
357             label faceI = myFaces[myFaceI];
359             if ((faceStatus[faceI] == ANYEDGE) || (faceStatus[faceI] >= 0))
360             {
361                 faceStatus[faceI] = NOEDGE;
362             }
363        }
364     }
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,
376     label edgeI
379     const edge& e = surf.edges()[edgeI];
380     label v1 = e.start();
381     label v2 = e.end();
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)
389     {
390         facesToBeCollapsed.insert(myFaces[myFaceI]);
391     }
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];
399     forAll(v1Faces, v1FaceI)
400     {
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);
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     }
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,
430     const label vertI
433     const labelList& myFaces = surf.pointFaces()[vertI];
435     forAll(myFaces, myFaceI)
436     {
437         label face1I = myFaces[myFaceI];
439         if (faceUsed.found(face1I))
440         {
441             return face1I;
442         }
443     }
444     return -1;
448 // Calculate new edge-face addressing resulting from edge collapse
449 void Foam::triSurfaceTools::getMergedEdges
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
458     const edge& e = surf.edges()[edgeI];
459     label v1 = e.start();
460     label v2 = e.end();
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)
469     {
470         if (!collapsedFaces.found(v2Faces[v2FaceI]))
471         {
472             v2FacesHash.insert(v2Faces[v2FaceI]);
473         }
474     }
477     forAll(v1Faces, v1FaceI)
478     {
479         label face1I = v1Faces[v1FaceI];
481         if (collapsedFaces.found(face1I))
482         {
483             continue;
484         }
486         //
487         // Check if face1I uses a vertex connected to a face using v2
488         //
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;
503         // Check vert1, vert2 for usage by v2Face.
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         }
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);
519             edgeToEdge.insert(edge1I, edge2I);
520             edgeToEdge.insert(edge2I, edge1I);
522             edgeToFace.insert(edge1I, face2I);
523             edgeToFace.insert(edge2I, face1I);
524         }
525     }
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,
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
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);
549     label D = -1;
551     label face2I = -1;
553     if (edgeToEdge.found(edgeI))
554     {
555         // Use merged addressing
556         label edge2I = edgeToEdge[edgeI];
557         face2I = edgeToFace[edgeI];
559         D = oppositeVertex(surf, face2I, edge2I);
560     }
561     else
562     {
563         // Use normal edge-face addressing
564         face2I = otherFace(surf, faceI, edgeI);
566         if ((face2I != -1) && !collapsedFaces.found(face2I))
567         {
568             D = oppositeVertex(surf, face2I, edgeI);
569         }
570     }
572     scalar cosAngle = 1;
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         {
618             FatalErrorIn("edgeCosAngle")
619                 << "face " << faceI << " does not use vertex "
620                 << v1 << " of collapsed edge" << abort(FatalError);
621         }
622     }
623     return cosAngle;
627 //- Calculate minimum (cos of) edge angle using addressing from collapsing
628 //  edge to v1 at pt.
629 Foam::scalar Foam::triSurfaceTools::collapseMinCosAngle
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
639     const labelList& v1Faces = surf.pointFaces()[v1];
641     scalar minCos = 1;
643     forAll(v1Faces, v1FaceI)
644     {
645         label faceI = v1Faces[v1FaceI];
647         if (collapsedFaces.found(faceI))
648         {
649             continue;
650         }
652         const labelList& myEdges = surf.faceEdges()[faceI];
654         forAll(myEdges, myEdgeI)
655         {
656             label edgeI = myEdges[myEdgeI];
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     }
677     return minCos;
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,
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
694     const labelList& v1Faces = surf.pointFaces()[v1];
696     forAll(v1Faces, v1FaceI)
697     {
698         label faceI = v1Faces[v1FaceI];
700         if (collapsedFaces.found(faceI))
701         {
702             continue;
703         }
705         const labelList& myEdges = surf.faceEdges()[faceI];
707         forAll(myEdges, myEdgeI)
708         {
709             label edgeI = myEdges[myEdgeI];
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     }
732     return false;
736 //// Return true if collapsing edgeI creates triangles on top of each other.
737 //// Is when the triangles neighbouring the collapsed one already share
738 // a vertex.
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)
759 //    {
760 //        const label faceI = collapsed[collapseI];
762 //        const labelList& myEdges = surf.faceEdges()[faceI];
764 //        Pout<< "collapsing faceI:" << faceI << " uses edges:" << myEdges
765 //            << endl;
767 //        forAll(myEdges, myEdgeI)
768 //        {
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))
775 //            {
776 //                // Get the other face
777 //                label neighbourFaceI = myFaces[0];
779 //                if (neighbourFaceI == faceI)
780 //                {
781 //                    neighbourFaceI = myFaces[1];
782 //                }
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 //    }
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.
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)
810 //    {
811 //        const labelList& faceIEdges = surf.faceEdges()[neighbourList[i]];
813 //        for (label j = i+1; j < neighbourList.size(); i++)
814 //        {
815 //            const labelList& faceJEdges = surf.faceEdges()[neighbourList[j]];
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;
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;
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()
842 Foam::surfaceLocation Foam::triSurfaceTools::cutEdge
844     const triSurface& s,
845     const label triI,
846     const label excludeEdgeI,
847     const label excludePointI,
849     const point& triPoint,
850     const plane& cutPlane,
851     const point& toPoint
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;
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     }
868     // Normalise and truncate
869     forAll(d, i)
870     {
871         d[i] /= norm;
873         if (mag(d[i]) < 1E-6)
874         {
875             d[i] = 0.0;
876         }
877     }
879     // Return information
880     surfaceLocation cut;
882     if (excludePointI != -1)
883     {
884         // Excluded point. Test only opposite edge.
886         label fp0 = findIndex(s.localFaces()[triI], excludePointI);
888         if (fp0 == -1)
889         {
890             FatalErrorIn("cutEdge(..)") << "excludePointI:" << excludePointI
891                 << " localF:" << s.localFaces()[triI] << abort(FatalError);
892         }
894         label fp1 = f.fcIndex(fp0);
895         label fp2 = f.fcIndex(fp1);
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
936         FixedList<surfaceLocation, 2> inters;
937         label interI = 0;
939         forAll(f, fp0)
940         {
941             label fp1 = f.fcIndex(fp0);
943             if (d[fp0] == 0)
944             {
945                 if (interI >= 2)
946                 {
947                     FatalErrorIn("cutEdge(..)")
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                 {
966                     FatalErrorIn("cutEdge(..)")
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         }
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;
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)
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()];
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()];
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()];
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()];
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     }
1161 // Start:
1162 //  - location
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);
1180     forAll(eFaces, i)
1181     {
1182         label triI = eFaces[i];
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.
1199                 surfaceLocation cutInfo = cutEdge
1200                 (
1201                     s,
1202                     triI,
1203                     excludeEdgeI,       // excludeEdgeI
1204                     excludePointI,      // excludePointI
1205                     start.rawPoint(),
1206                     cutPlane,
1207                     end.rawPoint()
1208                 );
1210                 // If crossing an edge we expect next edge to be cut.
1211                 if (excludeEdgeI != -1 && !cutInfo.hit())
1212                 {
1213                     FatalErrorIn("triSurfaceTools::visitFaces(..)")
1214                         << "Triangle:" << triI
1215                         << " excludeEdge:" << excludeEdgeI
1216                         << " point:" << start.rawPoint()
1217                         << " plane:" << cutPlane
1218                         << " . No intersection!" << abort(FatalError);
1219                 }
1221                 if (cutInfo.hit())
1222                 {
1223                     scalar distSqr = magSqr(cutInfo.rawPoint()-end.rawPoint());
1225                     if (distSqr < minDistSqr)
1226                     {
1227                         minDistSqr = distSqr;
1228                         nearest = cutInfo;
1229                         nearest.triangle() = triI;
1230                         nearest.setMiss();
1231                     }
1232                 }
1233             }
1234         }
1235     }
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     }
1244     return nearest;
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);
1260     forAll(pts, pointI)
1261     {
1262         const point& pt = pts[pointI];
1264         outFile<< "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
1265     }
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);
1280     label nVerts = 0;
1281     forAll(markedVerts, vertI)
1282     {
1283         if (markedVerts[vertI])
1284         {
1285             const point& pt = surf.localPoints()[vertI];
1287             outFile<< "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
1289             nVerts++;
1290         }
1291     }
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,
1304     const label edgeI,
1305     labelList& edgeTris
1308     const edge& e = surf.edges()[edgeI];
1309     const labelList& myFaces = surf.edgeFaces()[edgeI];
1311     label face1I = myFaces[0];
1312     label face2I = -1;
1313     if (myFaces.size() == 2)
1314     {
1315         face2I = myFaces[1];
1316     }
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());
1325     label nTris = 0;
1326     forAll(startFaces, startFaceI)
1327     {
1328         edgeTris[nTris++] = startFaces[startFaceI];
1329     }
1331     forAll(endFaces, endFaceI)
1332     {
1333         label faceI = endFaces[endFaceI];
1335         if ((faceI != face1I) && (faceI != face2I))
1336         {
1337             edgeTris[nTris++] = faceI;
1338         }
1339     }
1343 // Get all vertices connected to vertices of edge
1344 Foam::labelList Foam::triSurfaceTools::getVertexVertices
1346     const triSurface& surf,
1347     const edge& e
1350     const edgeList& edges = surf.edges();
1352     label v1 = e.start();
1353     label v2 = e.end();
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)
1361     {
1362         const edge& e = edges[v1Edges[v1EdgeI]];
1363         vertexNeighbours.insert(e.otherVertex(v1));
1364     }
1366     const labelList& v2Edges = surf.pointEdges()[v2];
1368     forAll(v2Edges, v2EdgeI)
1369     {
1370         const edge& e = edges[v2Edges[v2EdgeI]];
1372         label vertI = e.otherVertex(v2);
1374         vertexNeighbours.insert(vertI);
1375     }
1376     return vertexNeighbours.toc();
1380 //// Order vertices consistent with face
1381 //void Foam::triSurfaceTools::orderVertices
1383 //    const labelledTri& f,
1384 //    const label v1,
1385 //    const label v2,
1386 //    label& vA,
1387 //    label& vB
1390 //    // Order v1, v2 in anticlockwise order.
1391 //    bool reverse = false;
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 //    }
1415 //    if (reverse)
1416 //    {
1417 //        vA = v2;
1418 //        vB = v1;
1419 //    }
1420 //    else
1421 //    {
1422 //        vA = v1;
1423 //        vB = v2;
1424 //    }
1428 // Get the other face using edgeI
1429 Foam::label Foam::triSurfaceTools::otherFace
1431     const triSurface& surf,
1432     const label faceI,
1433     const label edgeI
1436     const labelList& myFaces = surf.edgeFaces()[edgeI];
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     }
1456 // Get the two edges on faceI counterclockwise after edgeI
1457 void Foam::triSurfaceTools::otherEdges
1459     const triSurface& surf,
1460     const label faceI,
1461     const label edgeI,
1462     label& e1,
1463     label& e2
1466     const labelList& eFaces = surf.faceEdges()[faceI];
1468     label i0 = findIndex(eFaces, edgeI);
1470     if (i0 == -1)
1471     {
1472         FatalErrorIn
1473         (
1474             "otherEdges"
1475             "(const triSurface&, const label, const label,"
1476             " label&, label&)"
1477         )   << "Edge " << surf.edges()[edgeI] << " not in face "
1478             << surf.localFaces()[faceI] << abort(FatalError);
1479     }
1481     label i1 = eFaces.fcIndex(i0);
1482     label i2 = eFaces.fcIndex(i1);
1484     e1 = eFaces[i1];
1485     e2 = eFaces[i2];
1489 // Get the two vertices on faceI counterclockwise vertI
1490 void Foam::triSurfaceTools::otherVertices
1492     const triSurface& surf,
1493     const label faceI,
1494     const label vertI,
1495     label& vert1I,
1496     label& vert2I
1499     const labelledTri& f = surf.localFaces()[faceI];
1501     if (vertI == f[0])
1502     {
1503         vert1I = f[1];
1504         vert2I = f[2];
1505     }
1506     else if (vertI == f[1])
1507     {
1508         vert1I = f[2];
1509         vert2I = f[0];
1510     }
1511     else if (vertI == f[2])
1512     {
1513         vert1I = f[0];
1514         vert2I = f[1];
1515     }
1516     else
1517     {
1518         FatalErrorIn
1519         (
1520             "otherVertices"
1521             "(const triSurface&, const label, const label,"
1522             " label&, label&)"
1523         )   << "Vertex " << vertI << " not in face " << f << abort(FatalError);
1524     }
1528 // Get edge opposite vertex
1529 Foam::label Foam::triSurfaceTools::oppositeEdge
1531     const triSurface& surf,
1532     const label faceI,
1533     const label vertI
1536     const labelList& myEdges = surf.faceEdges()[faceI];
1538     forAll(myEdges, myEdgeI)
1539     {
1540         label edgeI = myEdges[myEdgeI];
1542         const edge& e = surf.edges()[edgeI];
1544         if ((e.start() != vertI) && (e.end() != vertI))
1545         {
1546             return edgeI;
1547         }
1548     }
1550     FatalErrorIn
1551     (
1552         "oppositeEdge"
1553         "(const triSurface&, const label, const label)"
1554     )   << "Cannot find vertex " << vertI << " in edges of face " << faceI
1555         << abort(FatalError);
1557     return -1;
1561 // Get vertex opposite edge
1562 Foam::label Foam::triSurfaceTools::oppositeVertex
1564     const triSurface& surf,
1565     const label faceI,
1566     const label edgeI
1569     const triSurface::FaceType& f = surf.localFaces()[faceI];
1570     const edge& e = surf.edges()[edgeI];
1572     forAll(f, fp)
1573     {
1574         label vertI = f[fp];
1576         if (vertI != e.start() && vertI != e.end())
1577         {
1578             return vertI;
1579         }
1580     }
1582     FatalErrorIn("triSurfaceTools::oppositeVertex")
1583         << "Cannot find vertex opposite edge " << edgeI << " vertices " << e
1584         << " in face " << faceI << " vertices " << f << abort(FatalError);
1586     return -1;
1590 // Returns edge label connecting v1, v2
1591 Foam::label Foam::triSurfaceTools::getEdge
1593     const triSurface& surf,
1594     const label v1,
1595     const label v2
1598     const labelList& v1Edges = surf.pointEdges()[v1];
1600     forAll(v1Edges, v1EdgeI)
1601     {
1602         label edgeI = v1Edges[v1EdgeI];
1603         const edge& e = surf.edges()[edgeI];
1605         if ((e.start() == v2) || (e.end() == v2))
1606         {
1607             return edgeI;
1608         }
1609     }
1610     return -1;
1614 // Return index of triangle (or -1) using all three edges
1615 Foam::label Foam::triSurfaceTools::getTriangle
1617     const triSurface& surf,
1618     const label e0I,
1619     const label e1I,
1620     const label e2I
1623     if ((e0I == e1I) || (e0I == e2I) || (e1I == e2I))
1624     {
1625         FatalErrorIn
1626         (
1627             "getTriangle"
1628             "(const triSurface&, const label, const label,"
1629             " const label)"
1630         )   << "Duplicate edge labels : e0:" << e0I << " e1:" << e1I
1631             << " e2:" << e2I
1632             << abort(FatalError);
1633     }
1635     const labelList& eFaces = surf.edgeFaces()[e0I];
1637     forAll(eFaces, eFaceI)
1638     {
1639         label faceI = eFaces[eFaceI];
1641         const labelList& myEdges = surf.faceEdges()[faceI];
1643         if
1644         (
1645             (myEdges[0] == e1I)
1646          || (myEdges[1] == e1I)
1647          || (myEdges[2] == e1I)
1648         )
1649         {
1650             if
1651             (
1652                 (myEdges[0] == e2I)
1653              || (myEdges[1] == e2I)
1654              || (myEdges[2] == e2I)
1655             )
1656             {
1657                 return faceI;
1658             }
1659         }
1660     }
1661     return -1;
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)
1675     {
1676         const edge& e = surf.edges()[edgeI];
1678         edgeMids[edgeI] =
1679             0.5
1680           * (
1681                 surf.localPoints()[e.start()]
1682               + surf.localPoints()[e.end()]
1683             );
1684     }
1687     labelList faceStatus(surf.size(), ANYEDGE);
1689     //// Protect triangles which are on the border of different regions
1690     //forAll(edges, edgeI)
1691     //{
1692     //    const labelList& neighbours = edgeFaces[edgeI];
1693     //
1694     //    if ((neighbours.size() != 2) && (neighbours.size() != 1))
1695     //    {
1696     //        FatalErrorIn("collapseEdges")
1697     //            << abort(FatalError);
1698     //    }
1699     //
1700     //    if (neighbours.size() == 2)
1701     //    {
1702     //        if (surf[neighbours[0]].region() != surf[neighbours[1]].region())
1703     //        {
1704     //            // Neighbours on different regions. For now dont allow
1705     //            // any collapse.
1706     //            //Pout<< "protecting face " << neighbours[0]
1707     //            //    << ' ' << neighbours[1] << endl;
1708     //            faceStatus[neighbours[0]] = NOEDGE;
1709     //            faceStatus[neighbours[1]] = NOEDGE;
1710     //        }
1711     //    }
1712     //}
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)
1737     {
1738         pointMap[pointI] = pointI;
1739     }
1742     // Do actual 'collapsing' of edges
1744     forAll(collapseEdgeLabels, collapseEdgeI)
1745     {
1746         const label edgeI = collapseEdgeLabels[collapseEdgeI];
1748         if ((edgeI < 0) || (edgeI >= surf.nEdges()))
1749         {
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);
1755         }
1757         const labelList& neighbours = edgeFaces[edgeI];
1759         if (neighbours.size() == 2)
1760         {
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
1765             if
1766             (
1767                 ((stat0 == ANYEDGE) || (stat0 == edgeI))
1768              && ((stat1 == ANYEDGE) || (stat1 == edgeI))
1769             )
1770             {
1771                 const edge& e = edges[edgeI];
1773                 // Set up mapping to 'collapse' points of edge
1774                 if
1775                 (
1776                     (pointMap[e.start()] != e.start())
1777                  || (pointMap[e.end()] != e.end())
1778                 )
1779                 {
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);
1788                 }
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);
1800                 protectNeighbours
1801                 (
1802                     surf,
1803                     oppositeVertex(surf, neighbours[0], edgeI),
1804                     faceStatus
1805                 );
1806                 protectNeighbours
1807                 (
1808                     surf,
1809                     oppositeVertex(surf, neighbours[1], edgeI),
1810                     faceStatus
1811                 );
1813                 // Mark all collapsing faces
1814                 labelList collapseFaces =
1815                     getCollapsedFaces
1816                     (
1817                         surf,
1818                         edgeI
1819                     ).toc();
1821                 forAll(collapseFaces, collapseI)
1822                 {
1823                      faceStatus[collapseFaces[collapseI]] = COLLAPSED;
1824                 }
1825             }
1826         }
1827     }
1830     // Storage for new triangles
1831     List<labelledTri> newTris(surf.size());
1832     label newTriI = 0;
1834     const List<labelledTri>& localFaces = surf.localFaces();
1837     // Get only non-collapsed triangles and renumber vertex labels.
1838     forAll(localFaces, faceI)
1839     {
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]];
1846         if
1847         (
1848             (a != b) && (a != c) && (b != c)
1849          && (faceStatus[faceI] != COLLAPSED)
1850         )
1851         {
1852             // uncollapsed triangle
1853             newTris[newTriI++] = labelledTri(a, b, c, f.region());
1854         }
1855         else
1856         {
1857             //Pout<< "Collapsed triangle " << faceI
1858             //    << " vertices:" << f << endl;
1859         }
1860     }
1861     newTris.setSize(newTriI);
1865     // Pack faces
1867     triSurface tempSurf(newTris, surf.patches(), newPoints);
1869     return
1870         triSurface
1871         (
1872             tempSurf.localFaces(),
1873             tempSurf.patches(),
1874             tempSurf.localPoints()
1875         );
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)
1891     {
1892         calcRefineStatus(surf, refineFaces[refineFaceI], refineStatus);
1893     }
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();
1917     // Refine edges
1918     forAll(refineEdges, refineEdgeI)
1919     {
1920         label edgeI = refineEdges[refineEdgeI];
1922         const labelList& myFaces = surf.edgeFaces()[edgeI];
1924         bool neighbourIsRefined= false;
1926         forAll(myFaces, myFaceI)
1927         {
1928             if (refineStatus[myFaces[myFaceI]] != NONE)
1929             {
1930                 neighbourIsRefined =  true;
1931             }
1932         }
1934         // Only refine if none of the faces is refined
1935         if (!neighbourIsRefined)
1936         {
1937             // Refine edge
1938             const edge& e = surf.edges()[edgeI];
1940             point mid =
1941                 0.5
1942               * (
1943                     surf.localPoints()[e.start()]
1944                   + surf.localPoints()[e.end()]
1945                 );
1947             newPoints[newPointI] = mid;
1949             // Refine faces using edge
1950             forAll(myFaces, myFaceI)
1951             {
1952                 // Add faces to newFaces
1953                 greenRefine
1954                 (
1955                     surf,
1956                     myFaces[myFaceI],
1957                     edgeI,
1958                     newPointI,
1959                     newFaces
1960                 );
1962                 // Mark as refined
1963                 refineStatus[myFaces[myFaceI]] = GREEN;
1964             }
1966             newPointI++;
1967         }
1968     }
1970     // Add unrefined faces
1971     forAll(surf.localFaces(), faceI)
1972     {
1973         if (refineStatus[faceI] == NONE)
1974         {
1975             newFaces.append(surf.localFaces()[faceI]);
1976         }
1977     }
1979     newFaces.shrink();
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)
2001     {
2002         const edge& e = surf.edges()[edgeIndices[i]];
2004         scalar length =
2005             mag
2006             (
2007                 surf.localPoints()[e.end()]
2008               - surf.localPoints()[e.start()]
2009             );
2011         if (length < minLength)
2012         {
2013             minLength = length;
2014             minIndex = i;
2015         }
2016     }
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)
2031     {
2032         const edge& e = surf.edges()[edgeIndices[i]];
2034         scalar length =
2035             mag
2036             (
2037                 surf.localPoints()[e.end()]
2038               - surf.localPoints()[e.start()]
2039             );
2041         if (length > maxLength)
2042         {
2043             maxLength = length;
2044             maxIndex = i;
2045         }
2046     }
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
2063     (
2064         surf.localPoints(),
2065         mergeTol,
2066         false,
2067         pointMap,
2068         newPoints
2069     );
2071     if (hasMerged)
2072     {
2073         // Pack the triangles
2075         // Storage for new triangles
2076         List<labelledTri> newTriangles(surf.size());
2077         label newTriangleI = 0;
2079         forAll(surf, faceI)
2080         {
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))
2088             {
2089                 newTriangles[newTriangleI++] =
2090                     labelledTri(newA, newB, newC, f.region());
2091             }
2092         }
2093         newTriangles.setSize(newTriangleI);
2095         return triSurface
2096         (
2097             newTriangles,
2098             surf.patches(),
2099             newPoints
2100         );
2101     }
2102     else
2103     {
2104         return surf;
2105     }
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)
2125     {
2126         // Nearest to face
2127         return surf.faceNormals()[nearestFaceI];
2128     }
2129     else if (nearType == triPointRef::EDGE)
2130     {
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);
2139         forAll(eFaces, i)
2140         {
2141             edgeNormal += surf.faceNormals()[eFaces[i]];
2142         }
2143         return edgeNormal/(mag(edgeNormal) + VSMALL);
2144     }
2145     else
2146     {
2147         // Nearest to point
2148         const triSurface::FaceType& localF = surf.localFaces()[nearestFaceI];
2149         return surf.pointNormals()[localF[nearLabel]];
2150     }
2154 Foam::triSurfaceTools::sideType Foam::triSurfaceTools::edgeSide
2156     const triSurface& surf,
2157     const point& sample,
2158     const point& nearestPoint,
2159     const label edgeI
2162     const labelList& eFaces = surf.edgeFaces()[edgeI];
2164     if (eFaces.size() != 2)
2165     {
2166         // Surface not closed.
2167         return UNKNOWN;
2168     }
2169     else
2170     {
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)
2179         {
2180             return OUTSIDE;
2181         }
2182         else
2183         {
2184             return INSIDE;
2185         }
2186     }
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)
2209     {
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)
2221         // {
2222         //     c /= magSampleNearestVec*mag(surf.faceNormals()[nearestFaceI]);
2224         //     if (mag(c) < 0.99)
2225         //     {
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);
2243         //     }
2244         // }
2246         if (c > 0)
2247         {
2248             return OUTSIDE;
2249         }
2250         else
2251         {
2252             return INSIDE;
2253         }
2254     }
2255     else if (nearType == triPointRef::EDGE)
2256     {
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];
2263         // if (debug)
2264         // {
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]]);
2270         //    if
2271         //    (
2272         //        meshEdge
2273         //     != edge(f[nearLabel], f[f.fcIndex(nearLabel)])
2274         //    )
2275         //    {
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);
2282         //    }
2283         // }
2285         return edgeSide(surf, sample, nearestPoint, edgeI);
2286     }
2287     else
2288     {
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;
2306         forAll(pEdges, i)
2307         {
2308             label edgeI = pEdges[i];
2310             const edge& e = edges[edgeI];
2312             label otherPointI = e.otherVertex(nearPointI);
2314             // Get edge normal.
2315             vector eVec(localPoints[otherPointI] - base);
2316             scalar magEVec = mag(eVec);
2318             if (magEVec > VSMALL)
2319             {
2320                 eVec /= magEVec;
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)
2328                 {
2329                     minDistSqr = distSqr;
2330                     minEdgeI = edgeI;
2331                 }
2332             }
2333         }
2335         if (minEdgeI == -1)
2336         {
2337             FatalErrorIn("treeDataTriSurface::getSide")
2338                 << "Problem: did not find edge closer than " << minDistSqr
2339                 << abort(FatalError);
2340         }
2342         return edgeSide(surf, sample, nearestPoint, minEdgeI);
2343     }
2347 // triangulation of boundaryMesh
2348 Foam::triSurface Foam::triSurfaceTools::triangulate
2350     const polyBoundaryMesh& bMesh,
2351     const labelHashSet& includePatches,
2352     const bool verbose
2355     const polyMesh& mesh = bMesh.mesh();
2357     // Storage for surfaceMesh. Size estimate.
2358     DynamicList<labelledTri> triangles
2359     (
2360         mesh.nFaces() - mesh.nInternalFaces()
2361     );
2363     label newPatchI = 0;
2365     forAllConstIter(labelHashSet, includePatches, iter)
2366     {
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)
2374         {
2375             const face& f = patch[patchFaceI];
2377             faceList triFaces(f.nTriangles(points));
2379             label nTri = 0;
2381             f.triangles(points, nTri, triFaces);
2383             forAll(triFaces, triFaceI)
2384             {
2385                 const face& f = triFaces[triFaceI];
2387                 triangles.append(labelledTri(f[0], f[1], f[2], newPatchI));
2389                 nTriTotal++;
2390             }
2391         }
2393         if (verbose)
2394         {
2395             Pout<< patch.name() << " : generated " << nTriTotal
2396                 << " triangles from " << patch.size() << " faces with"
2397                 << " new patchid " << newPatchI << endl;
2398         }
2400         newPatchI++;
2401     }
2402     triangles.shrink();
2404     // Create globally numbered tri surface
2405     triSurface rawSurface(triangles, mesh.points());
2407     // Create locally numbered tri surface
2408     triSurface surface
2409     (
2410         rawSurface.localFaces(),
2411         rawSurface.localPoints()
2412     );
2414     // Add patch names to surface
2415     surface.patches().setSize(newPatchI);
2417     newPatchI = 0;
2419     forAllConstIter(labelHashSet, includePatches, iter)
2420     {
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();
2427         newPatchI++;
2428     }
2430     return surface;
2434 // triangulation of boundaryMesh
2435 Foam::triSurface Foam::triSurfaceTools::triangulateFaceCentre
2437     const polyBoundaryMesh& bMesh,
2438     const labelHashSet& includePatches,
2439     const bool verbose
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)
2453     {
2454         newPoints[newPointI++] = points[pointI];
2455     }
2456     forAll(faceCentres, faceI)
2457     {
2458         newPoints[newPointI++] = faceCentres[faceI];
2459     }
2462     // Count number of faces.
2463     DynamicList<labelledTri> triangles
2464     (
2465         mesh.nFaces() - mesh.nInternalFaces()
2466     );
2468     label newPatchI = 0;
2470     forAllConstIter(labelHashSet, includePatches, iter)
2471     {
2472         const label patchI = iter.key();
2473         const polyPatch& patch = bMesh[patchI];
2475         label nTriTotal = 0;
2477         forAll(patch, patchFaceI)
2478         {
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();
2485             forAll(f, fp)
2486             {
2487                 label fp1 = f.fcIndex(fp);
2489                 triangles.append(labelledTri(f[fp], f[fp1], fc, newPatchI));
2491                 nTriTotal++;
2492             }
2493         }
2495         if (verbose)
2496         {
2497             Pout<< patch.name() << " : generated " << nTriTotal
2498                 << " triangles from " << patch.size() << " faces with"
2499                 << " new patchid " << newPatchI << endl;
2500         }
2502         newPatchI++;
2503     }
2504     triangles.shrink();
2507     // Create globally numbered tri surface
2508     triSurface rawSurface(triangles, newPoints);
2510     // Create locally numbered tri surface
2511     triSurface surface
2512     (
2513         rawSurface.localFaces(),
2514         rawSurface.localPoints()
2515     );
2517     // Add patch names to surface
2518     surface.patches().setSize(newPatchI);
2520     newPatchI = 0;
2522     forAllConstIter(labelHashSet, includePatches, iter)
2523     {
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();
2530         newPatchI++;
2531     }
2533     return surface;
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());
2542     label doubleI = 0;
2543     forAll(pts, i)
2544     {
2545         geompackVertices[doubleI++] = pts[i][0];
2546         geompackVertices[doubleI++] = pts[i][1];
2547     }
2549     // Storage for triangles
2550     int m2 = 3;
2551     List<int> triangle_node(m2*3*pts.size());
2552     List<int> triangle_neighbor(m2*3*pts.size());
2554     // Triangulate
2555     int nTris = 0;
2556     int err = dtris2
2557     (
2558         pts.size(),
2559         geompackVertices.begin(),
2560         &nTris,
2561         triangle_node.begin(),
2562         triangle_neighbor.begin()
2563     );
2565     if (err != 0)
2566     {
2567         FatalErrorIn("triSurfaceTools::delaunay2D(const List<vector2D>&)")
2568             << "Failed dtris2 with vertices:" << pts.size()
2569             << abort(FatalError);
2570     }
2572     // Trim
2573     triangle_node.setSize(3*nTris);
2574     triangle_neighbor.setSize(3*nTris);
2576     // Convert to triSurface.
2577     List<labelledTri> faces(nTris);
2579     forAll(faces, i)
2580     {
2581         faces[i] = labelledTri
2582         (
2583             triangle_node[3*i]-1,
2584             triangle_node[3*i+1]-1,
2585             triangle_node[3*i+2]-1,
2586             0
2587         );
2588     }
2590     pointField points(pts.size());
2591     forAll(pts, i)
2592     {
2593         points[i][0] = pts[i][0];
2594         points[i][1] = pts[i][1];
2595         points[i][2] = 0.0;
2596     }
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)
2621 //    Triangulation T;
2623 //    // Insert 2D vertices; building triangulation
2624 //    forAll(pts, i)
2625 //    {
2626 //        const point& pt = pts[i];
2628 //        T.insert(Point(pt[0], pt[1]));
2629 //    }
2632 //    // Number vertices
2633 //    // ~~~~~~~~~~~~~~~
2635 //    label vertI = 0;
2636 //    for
2637 //    (
2638 //        Vertex_iterator it = T.vertices_begin();
2639 //        it != T.vertices_end();
2640 //        ++it
2641 //    )
2642 //    {
2643 //        it->info() = vertI++;
2644 //    }
2646 //    // Extract faces
2647 //    // ~~~~~~~~~~~~~
2649 //    List<labelledTri> faces(T.number_of_faces());
2651 //    label faceI = 0;
2653 //    for
2654 //    (
2655 //        Finite_faces_iterator fc = T.finite_faces_begin();
2656 //        fc != T.finite_faces_end();
2657 //        ++fc
2658 //    )
2659 //    {
2660 //        faces[faceI++] = labelledTri
2661 //        (
2662 //            fc->vertex(0)->info(),
2663 //            f[1] = fc->vertex(1)->info(),
2664 //            f[2] = fc->vertex(2)->info()
2665 //        );
2666 //    }
2668 //    pointField points(pts.size());
2669 //    forAll(pts, i)
2670 //    {
2671 //        points[i][0] = pts[i][0];
2672 //        points[i][1] = pts[i][1];
2673 //        points[i][2] = 0.0;
2674 //    }
2676 //    return triSurface(faces, points);
2680 void Foam::triSurfaceTools::calcInterpolationWeights
2682     const triPointRef& tri,
2683     const point& p,
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++)
2699     {
2700         normal[i] = triangleFaceNormal ^ edge[i];
2701         normal[i] /= mag(normal[i]) + VSMALL;
2702     }
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)
2726     {
2727         const point& samplePt = samplePts[i];
2730         FixedList<label, 3>& verts = allVerts[i];
2731         FixedList<scalar, 3>& weights = allWeights[i];
2733         scalar minDistance = GREAT;
2735         forAll(s, faceI)
2736         {
2737             const labelledTri& f = s[faceI];
2739             triPointRef tri(f.tri(points));
2741             label nearType, nearLabel;
2743             pointHit nearest = tri.nearestPointClassify
2744             (
2745                 samplePt,
2746                 nearType,
2747                 nearLabel
2748             );
2750             if (nearest.hit())
2751             {
2752                 // samplePt inside triangle
2753                 verts[0] = f[0];
2754                 verts[1] = f[1];
2755                 verts[2] = f[2];
2757                 calcInterpolationWeights(tri, nearest.rawPoint(), weights);
2759                 //Pout<< "calcScalingFactors : samplePt:" << samplePt
2760                 //    << " inside triangle:" << faceI
2761                 //    << " verts:" << verts
2762                 //    << " weights:" << weights
2763                 //    << endl;
2765                 break;
2766             }
2767             else if (nearest.distance() < minDistance)
2768             {
2769                 minDistance = nearest.distance();
2771                 // Outside triangle. Store nearest.
2773                 if (nearType == triPointRef::POINT)
2774                 {
2775                     verts[0] = f[nearLabel];
2776                     weights[0] = 1;
2777                     verts[1] = -1;
2778                     weights[1] = -GREAT;
2779                     verts[2] = -1;
2780                     weights[2] = -GREAT;
2782                     //Pout<< "calcScalingFactors : samplePt:" << samplePt
2783                     //    << " distance:" << nearest.distance()
2784                     //    << " from point:" << points[f[nearLabel]]
2785                     //    << endl;
2786                 }
2787                 else if (nearType == triPointRef::EDGE)
2788                 {
2789                     verts[0] = f[nearLabel];
2790                     verts[1] = f[f.fcIndex(nearLabel)];
2791                     verts[2] = -1;
2793                     const point& p0 = points[verts[0]];
2794                     const point& p1 = points[verts[1]];
2796                     scalar s = min
2797                     (
2798                         1,
2799                         max
2800                         (
2801                             0,
2802                             mag(nearest.rawPoint() - p0)/mag(p1 - p0)
2803                         )
2804                     );
2806                     // Interpolate
2807                     weights[0] = 1 - s;
2808                     weights[1] = s;
2809                     weights[2] = -GREAT;
2811                     //Pout<< "calcScalingFactors : samplePt:" << samplePt
2812                     //    << " distance:" << nearest.distance()
2813                     //    << " from edge:" << p0 << p1 << " s:" << s
2814                     //    << endl;
2815                 }
2816                 else
2817                 {
2818                     // triangle. Can only happen because of truncation errors.
2819                     verts[0] = f[0];
2820                     verts[1] = f[1];
2821                     verts[2] = f[2];
2823                     calcInterpolationWeights(tri, nearest.rawPoint(), weights);
2825                     break;
2826                 }
2827             }
2828         }
2829     }
2833 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2834 // Tracking:
2837 // Test point on surface to see if is on face,edge or point.
2838 Foam::surfaceLocation Foam::triSurfaceTools::classify
2840     const triSurface& s,
2841     const label triI,
2842     const point& trianglePoint
2845     surfaceLocation nearest;
2847     // Nearest point could be on point or edge. Retest.
2848     label index, elemType;
2849     //bool inside =
2850     triPointRef(s[triI].tri(s.points())).classify
2851     (
2852         trianglePoint,
2853         elemType,
2854         index
2855     );
2857     nearest.setPoint(trianglePoint);
2859     if (elemType == triPointRef::NONE)
2860     {
2861         nearest.setHit();
2862         nearest.setIndex(triI);
2863         nearest.elementType() = triPointRef::NONE;
2864     }
2865     else if (elemType == triPointRef::EDGE)
2866     {
2867         nearest.setMiss();
2868         nearest.setIndex(s.faceEdges()[triI][index]);
2869         nearest.elementType() = triPointRef::EDGE;
2870     }
2871     else //if (elemType == triPointRef::POINT)
2872     {
2873         nearest.setMiss();
2874         nearest.setIndex(s.localFaces()[triI][index]);
2875         nearest.elementType() = triPointRef::POINT;
2876     }
2878     return nearest;
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;
2892     nearest.setMiss();
2894     // See if in same triangle as endpoint. If so snap.
2895     snapToEnd(s, end, nearest);
2897     if (!nearest.hit())
2898     {
2899         // Not yet at end point
2901         if (start.elementType() == triPointRef::NONE)
2902         {
2903             // Start point is inside triangle. Trivial cases already handled
2904             // above.
2906             // end point is on edge or point so cross currrent triangle to
2907             // see which edge is cut.
2909             nearest = cutEdge
2910             (
2911                 s,
2912                 start.index(),          // triangle
2913                 -1,                     // excludeEdge
2914                 -1,                     // excludePoint
2915                 start.rawPoint(),
2916                 cutPlane,
2917                 end.rawPoint()
2918             );
2919             nearest.elementType() = triPointRef::EDGE;
2920             nearest.triangle() = start.index();
2921             nearest.setMiss();
2922         }
2923         else if (start.elementType() == triPointRef::EDGE)
2924         {
2925             // Pick connected triangle that is most in direction.
2926             const labelList& eFaces = s.edgeFaces()[start.index()];
2928             nearest = visitFaces
2929             (
2930                 s,
2931                 eFaces,
2932                 start,
2933                 start.index(),      // excludeEdgeI
2934                 -1,                 // excludePointI
2935                 end,
2936                 cutPlane
2937             );
2938         }
2939         else    // start.elementType() == triPointRef::POINT
2940         {
2941             const labelList& pFaces = s.pointFaces()[start.index()];
2943             nearest = visitFaces
2944             (
2945                 s,
2946                 pFaces,
2947                 start,
2948                 -1,                 // excludeEdgeI
2949                 start.index(),      // excludePointI
2950                 end,
2951                 cutPlane
2952             );
2953         }
2954         snapToEnd(s, end, nearest);
2955     }
2956     return 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");
2969     //label vertI = 0;
2970     //meshTools::writeOBJ(str, hitInfo.rawPoint());
2971     //vertI++;
2973     // Track across surface.
2974     while (true)
2975     {
2976         //Pout<< "Tracking from:" << nl
2977         //    << "    " << hitInfo.info()
2978         //    << endl;
2980         hitInfo = trackToEdge
2981         (
2982             s,
2983             hitInfo,
2984             endInfo,
2985             cutPlane
2986         );
2988         //meshTools::writeOBJ(str, hitInfo.rawPoint());
2989         //vertI++;
2990         //str<< "l " << vertI-1 << ' ' << vertI << nl;
2992         //Pout<< "Tracked to:" << nl
2993         //    << "    " << hitInfo.info() << endl;
2995         if (hitInfo.hit() || hitInfo.triangle() == -1)
2996         {
2997             break;
2998         }
2999     }
3003 // ************************************************************************* //