triSurfaceTools.C 77 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
OpenFOAM bot's avatar
OpenFOAM bot committed
5
    \\  /    A nd           | www.openfoam.com
OpenFOAM bot's avatar
OpenFOAM bot committed
6 7
     \\/     M anipulation  |
-------------------------------------------------------------------------------
OpenFOAM bot's avatar
OpenFOAM bot committed
8
    Copyright (C) 2011-2016 OpenFOAM Foundation
9
    Copyright (C) 2016-2020 OpenCFD Ltd.
10 11 12 13
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

14 15 16 17
    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.
18 19 20 21 22 23 24

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
25
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
26 27 28 29 30 31

\*---------------------------------------------------------------------------*/

#include "triSurfaceTools.H"

#include "triSurface.H"
32
#include "MeshedSurface.H"
33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54
#include "OFstream.H"
#include "mergePoints.H"
#include "polyMesh.H"
#include "plane.H"
#include "geompack.H"

// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //

const Foam::label Foam::triSurfaceTools::ANYEDGE = -1;
const Foam::label Foam::triSurfaceTools::NOEDGE = -2;
const Foam::label Foam::triSurfaceTools::COLLAPSED = -3;

// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //

/*
    Refine by splitting all three edges of triangle ('red' refinement).
    Neighbouring triangles (which are not marked for refinement get split
    in half ('green') refinement. (R. Verfuerth, "A review of a posteriori
    error estimation and adaptive mesh refinement techniques",
    Wiley-Teubner, 1996)
*/

55
// Facei gets refined ('red'). Propagate edge refinement.
56 57 58
void Foam::triSurfaceTools::calcRefineStatus
(
    const triSurface& surf,
59
    const label facei,
60 61 62
    List<refineType>& refine
)
{
63
    if (refine[facei] == RED)
64 65 66 67 68 69
    {
        // Already marked for refinement. Do nothing.
    }
    else
    {
        // Not marked or marked for 'green' refinement. Refine.
70
        refine[facei] = RED;
71

72
        const labelList& myNeighbours = surf.faceFaces()[facei];
73

74
        for (const label neighbourFacei : myNeighbours)
75
        {
76
            if (refine[neighbourFacei] == GREEN)
77 78
            {
                // Change to red refinement and propagate
79
                calcRefineStatus(surf, neighbourFacei, refine);
80
            }
81
            else if (refine[neighbourFacei] == NONE)
82
            {
83
                refine[neighbourFacei] = GREEN;
84 85 86 87 88 89
            }
        }
    }
}


90
// Split facei along edgeI at position newPointi
91 92 93
void Foam::triSurfaceTools::greenRefine
(
    const triSurface& surf,
94
    const label facei,
95
    const label edgeI,
96
    const label newPointi,
97 98 99
    DynamicList<labelledTri>& newFaces
)
{
100
    const labelledTri& f = surf.localFaces()[facei];
101 102 103 104
    const edge& e = surf.edges()[edgeI];

    // Find index of edge in face.

105
    label fp0 = f.find(e[0]);
106 107 108 109 110 111 112 113 114 115 116
    label fp1 = f.fcIndex(fp0);
    label fp2 = f.fcIndex(fp1);

    if (f[fp1] == e[1])
    {
        // Edge oriented like face
        newFaces.append
        (
            labelledTri
            (
                f[fp0],
117
                newPointi,
118 119 120 121 122 123 124 125
                f[fp2],
                f.region()
            )
        );
        newFaces.append
        (
            labelledTri
            (
126
                newPointi,
127 128 129 130 131 132 133 134 135 136 137 138 139
                f[fp1],
                f[fp2],
                f.region()
            )
        );
    }
    else
    {
        newFaces.append
        (
            labelledTri
            (
                f[fp2],
140
                newPointi,
141 142 143 144 145 146 147 148
                f[fp1],
                f.region()
            )
        );
        newFaces.append
        (
            labelledTri
            (
149
                newPointi,
150 151 152 153 154 155 156 157 158
                f[fp0],
                f[fp1],
                f.region()
            )
        );
    }
}


159
// Refine all triangles marked for refinement.
160 161 162 163 164 165 166 167 168
Foam::triSurface Foam::triSurfaceTools::doRefine
(
    const triSurface& surf,
    const List<refineType>& refineStatus
)
{
    // Storage for new points. (start after old points)
    label newVertI = surf.nPoints();

169 170 171
    DynamicList<point> newPoints(newVertI);
    newPoints.append(surf.localPoints());

172 173 174 175 176 177 178
    // Storage for new faces
    DynamicList<labelledTri> newFaces(surf.size());


    // Point index for midpoint on edge
    labelList edgeMid(surf.nEdges(), -1);

179
    forAll(refineStatus, facei)
180
    {
181
        if (refineStatus[facei] == RED)
182 183
        {
            // Create new vertices on all edges to be refined.
184
            const labelList& fEdges = surf.faceEdges()[facei];
185

186
            for (const label edgei : fEdges)
187
            {
188
                if (edgeMid[edgei] == -1)
189
                {
190
                    const edge& e = surf.edges()[edgei];
191 192

                    // Create new point on mid of edge
193 194
                    newPoints.append(e.centre(surf.localPoints()));
                    edgeMid[edgei] = newVertI++;
195 196 197 198
                }
            }

            // Now we have new mid edge vertices for all edges on face
Andrew Heather's avatar
Andrew Heather committed
199
            // so create triangles for RED refinement.
200 201 202 203 204 205 206 207 208 209 210

            const edgeList& edges = surf.edges();

            // Corner triangles
            newFaces.append
            (
                labelledTri
                (
                    edgeMid[fEdges[0]],
                    edges[fEdges[0]].commonVertex(edges[fEdges[1]]),
                    edgeMid[fEdges[1]],
211
                    surf[facei].region()
212 213 214 215 216 217 218 219 220 221
                )
            );

            newFaces.append
            (
                labelledTri
                (
                    edgeMid[fEdges[1]],
                    edges[fEdges[1]].commonVertex(edges[fEdges[2]]),
                    edgeMid[fEdges[2]],
222
                    surf[facei].region()
223 224 225 226 227 228 229 230 231 232
                )
            );

            newFaces.append
            (
                labelledTri
                (
                    edgeMid[fEdges[2]],
                    edges[fEdges[2]].commonVertex(edges[fEdges[0]]),
                    edgeMid[fEdges[0]],
233
                    surf[facei].region()
234 235 236 237 238 239 240 241 242 243 244
                )
            );

            // Inner triangle
            newFaces.append
            (
                labelledTri
                (
                    edgeMid[fEdges[0]],
                    edgeMid[fEdges[1]],
                    edgeMid[fEdges[2]],
245
                    surf[facei].region()
246 247 248 249 250
                )
            );


            // Create triangles for GREEN refinement.
251
            for (const label edgei : fEdges)
252
            {
253
                label otherFacei = otherFace(surf, facei, edgei);
254

255
                if ((otherFacei != -1) && (refineStatus[otherFacei] == GREEN))
256 257 258 259
                {
                    greenRefine
                    (
                        surf,
260
                        otherFacei,
261 262
                        edgei,
                        edgeMid[edgei],
263 264 265 266 267 268 269 270
                        newFaces
                    );
                }
            }
        }
    }

    // Copy unmarked triangles since keep original vertex numbering.
271
    forAll(refineStatus, facei)
272
    {
273
        if (refineStatus[facei] == NONE)
274
        {
275
            newFaces.append(surf.localFaces()[facei]);
276 277 278 279 280 281 282 283 284 285 286 287
        }
    }

    newFaces.shrink();
    newPoints.shrink();


    // Transfer DynamicLists to straight ones.
    pointField allPoints;
    allPoints.transfer(newPoints);
    newPoints.clear();

mattijs's avatar
mattijs committed
288
    return triSurface(newFaces, surf.patches(), allPoints, true);
289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305
}


// Angle between two neighbouring triangles,
// angle between shared-edge end points and left and right face end points
Foam::scalar Foam::triSurfaceTools::faceCosAngle
(
    const point& pStart,
    const point& pEnd,
    const point& pLeft,
    const point& pRight
)
{
    const vector common(pEnd - pStart);
    const vector base0(pLeft - pStart);
    const vector base1(pRight - pStart);

306 307
    const vector n0 = normalised(common ^ base0);
    const vector n1 = normalised(base1 ^ common);
308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323

    return n0 & n1;
}


// Protect edges around vertex from collapsing.
// Moving vertI to a different position
// - affects obviously the cost of the faces using it
// - but also their neighbours since the edge between the faces changes
void Foam::triSurfaceTools::protectNeighbours
(
    const triSurface& surf,
    const label vertI,
    labelList& faceStatus
)
{
324
//    for (const label facei : surf.pointFaces()[vertI])
325
//    {
326
//        if ((faceStatus[facei] == ANYEDGE) || (faceStatus[facei] >= 0))
327
//        {
328
//            faceStatus[facei] = NOEDGE;
329 330 331
//        }
//    }

332
    for (const label edgei : surf.pointEdges()[vertI])
333
    {
334
        for (const label facei : surf.edgeFaces()[edgei])
335
        {
336
            if ((faceStatus[facei] == ANYEDGE) || (faceStatus[facei] >= 0))
337
            {
338
                faceStatus[facei] = NOEDGE;
339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356
            }
       }
    }
}


//
// Edge collapse helper functions
//

// Get all faces that will get collapsed if edgeI collapses.
Foam::labelHashSet Foam::triSurfaceTools::getCollapsedFaces
(
    const triSurface& surf,
    label edgeI
)
{
    const edge& e = surf.edges()[edgeI];
357 358
    const label v1 = e.start();
    const label v2 = e.end();
359 360 361 362 363

    // Faces using edge will certainly get collapsed.
    const labelList& myFaces = surf.edgeFaces()[edgeI];

    labelHashSet facesToBeCollapsed(2*myFaces.size());
364
    facesToBeCollapsed.insert(myFaces);
365

366 367
    // From faces using v1 check if they share an edge with faces
    // using v2.
368
    //  - share edge: are part of 'splay' tree and will collapse if edge
369 370 371
    //    collapses
    const labelList& v1Faces = surf.pointFaces()[v1];

372
    for (const label face1I : v1Faces)
373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403
    {
        label otherEdgeI = oppositeEdge(surf, face1I, v1);

        // Step across edge to other face
        label face2I = otherFace(surf, face1I, otherEdgeI);

        if (face2I != -1)
        {
            // found face on other side of edge. Now check if uses v2.
            if (oppositeVertex(surf, face2I, otherEdgeI) == v2)
            {
                // triangles face1I and face2I form splay tree and will
                // collapse.
                facesToBeCollapsed.insert(face1I);
                facesToBeCollapsed.insert(face2I);
            }
        }
    }

    return facesToBeCollapsed;
}


// Return value of faceUsed for faces using vertI
Foam::label Foam::triSurfaceTools::vertexUsesFace
(
    const triSurface& surf,
    const labelHashSet& faceUsed,
    const label vertI
)
{
404
    for (const label face1I : surf.pointFaces()[vertI])
405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420
    {
        if (faceUsed.found(face1I))
        {
            return face1I;
        }
    }
    return -1;
}


// Calculate new edge-face addressing resulting from edge collapse
void Foam::triSurfaceTools::getMergedEdges
(
    const triSurface& surf,
    const label edgeI,
    const labelHashSet& collapsedFaces,
421 422
    Map<label>& edgeToEdge,
    Map<label>& edgeToFace
423 424 425
)
{
    const edge& e = surf.edges()[edgeI];
426 427
    const label v1 = e.start();
    const label v2 = e.end();
428 429 430 431 432 433 434

    const labelList& v1Faces = surf.pointFaces()[v1];
    const labelList& v2Faces = surf.pointFaces()[v2];

    // Mark all (non collapsed) faces using v2
    labelHashSet v2FacesHash(v2Faces.size());

435
    for (const label facei : v2Faces)
436
    {
437
        if (!collapsedFaces.found(facei))
438
        {
439
            v2FacesHash.insert(facei);
440 441 442 443
        }
    }


444
    for (const label face1I: v1Faces)
445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493
    {
        if (collapsedFaces.found(face1I))
        {
            continue;
        }

        //
        // Check if face1I uses a vertex connected to a face using v2
        //

        label vert1I = -1;
        label vert2I = -1;
        otherVertices
        (
            surf,
            face1I,
            v1,
            vert1I,
            vert2I
        );
        //Pout<< "Face:" << surf.localFaces()[face1I] << " other vertices:"
        //    << vert1I << ' ' << vert2I << endl;

        // Check vert1, vert2 for usage by v2Face.

        label commonVert = vert1I;
        label face2I = vertexUsesFace(surf, v2FacesHash, commonVert);
        if (face2I == -1)
        {
            commonVert = vert2I;
            face2I = vertexUsesFace(surf, v2FacesHash, commonVert);
        }

        if (face2I != -1)
        {
            // Found one: commonVert is used by both face1I and face2I
            label edge1I = getEdge(surf, v1, commonVert);
            label edge2I = getEdge(surf, v2, commonVert);

            edgeToEdge.insert(edge1I, edge2I);
            edgeToEdge.insert(edge2I, edge1I);

            edgeToFace.insert(edge1I, face2I);
            edgeToFace.insert(edge2I, face1I);
        }
    }
}


494
// Calculates (cos of) angle across edgeI of facei,
495 496 497 498 499 500 501
// taking into account updated addressing (resulting from edge collapse)
Foam::scalar Foam::triSurfaceTools::edgeCosAngle
(
    const triSurface& surf,
    const label v1,
    const point& pt,
    const labelHashSet& collapsedFaces,
502 503
    const Map<label>& edgeToEdge,
    const Map<label>& edgeToFace,
504
    const label facei,
505 506 507 508 509 510 511
    const label edgeI
)
{
    const pointField& localPoints = surf.localPoints();

    label A = surf.edges()[edgeI].start();
    label B = surf.edges()[edgeI].end();
512
    label C = oppositeVertex(surf, facei, edgeI);
513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528

    label D = -1;

    label face2I = -1;

    if (edgeToEdge.found(edgeI))
    {
        // Use merged addressing
        label edge2I = edgeToEdge[edgeI];
        face2I = edgeToFace[edgeI];

        D = oppositeVertex(surf, face2I, edge2I);
    }
    else
    {
        // Use normal edge-face addressing
529
        face2I = otherFace(surf, facei, edgeI);
530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582

        if ((face2I != -1) && !collapsedFaces.found(face2I))
        {
            D = oppositeVertex(surf, face2I, edgeI);
        }
    }

    scalar cosAngle = 1;

    if (D != -1)
    {
        if (A == v1)
        {
            cosAngle = faceCosAngle
            (
                pt,
                localPoints[B],
                localPoints[C],
                localPoints[D]
            );
        }
        else if (B == v1)
        {
            cosAngle = faceCosAngle
            (
                localPoints[A],
                pt,
                localPoints[C],
                localPoints[D]
            );
        }
        else if (C == v1)
        {
            cosAngle = faceCosAngle
            (
                localPoints[A],
                localPoints[B],
                pt,
                localPoints[D]
            );
        }
        else if (D == v1)
        {
            cosAngle = faceCosAngle
            (
                localPoints[A],
                localPoints[B],
                localPoints[C],
                pt
            );
        }
        else
        {
583
            FatalErrorInFunction
584
                << "face " << facei << " does not use vertex "
585 586 587 588 589 590 591 592 593 594 595 596 597
                << v1 << " of collapsed edge" << abort(FatalError);
        }
    }
    return cosAngle;
}


Foam::scalar Foam::triSurfaceTools::collapseMinCosAngle
(
    const triSurface& surf,
    const label v1,
    const point& pt,
    const labelHashSet& collapsedFaces,
598 599
    const Map<label>& edgeToEdge,
    const Map<label>& edgeToFace
600 601 602 603 604 605
)
{
    const labelList& v1Faces = surf.pointFaces()[v1];

    scalar minCos = 1;

606
    for (const label facei : v1Faces)
607
    {
608
        if (collapsedFaces.found(facei))
609 610 611 612
        {
            continue;
        }

613
        for (const label edgeI : surf.faceEdges()[facei])
614 615 616 617 618 619 620 621 622 623 624 625 626
        {
            minCos =
                min
                (
                    minCos,
                    edgeCosAngle
                    (
                        surf,
                        v1,
                        pt,
                        collapsedFaces,
                        edgeToEdge,
                        edgeToFace,
627
                        facei,
628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645
                        edgeI
                    )
                );
        }
    }

    return minCos;
}


// Calculate max dihedral angle after collapsing edge to v1 (at pt).
// Note that all edges of all faces using v1 are affected.
bool Foam::triSurfaceTools::collapseCreatesFold
(
    const triSurface& surf,
    const label v1,
    const point& pt,
    const labelHashSet& collapsedFaces,
646 647
    const Map<label>& edgeToEdge,
    const Map<label>& edgeToFace,
648 649 650 651 652
    const scalar minCos
)
{
    const labelList& v1Faces = surf.pointFaces()[v1];

653
    for (const label facei : v1Faces)
654
    {
655
        if (collapsedFaces.found(facei))
656 657 658 659
        {
            continue;
        }

660
        const labelList& myEdges = surf.faceEdges()[facei];
661

662
        for (const label edgeI : myEdges)
663 664 665 666 667 668 669 670 671 672 673
        {
            if
            (
                edgeCosAngle
                (
                    surf,
                    v1,
                    pt,
                    collapsedFaces,
                    edgeToEdge,
                    edgeToFace,
674
                    facei,
675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706
                    edgeI
                )
              < minCos
            )
            {
                return true;
            }
        }
    }

    return false;
}


//// Return true if collapsing edgeI creates triangles on top of each other.
//// Is when the triangles neighbouring the collapsed one already share
// a vertex.
//bool Foam::triSurfaceTools::collapseCreatesDuplicates
//(
//    const triSurface& surf,
//    const label edgeI,
//    const labelHashSet& collapsedFaces
//)
//{
//Pout<< "duplicate : edgeI:" << edgeI << surf.edges()[edgeI]
//    << " collapsedFaces:" << collapsedFaces.toc() << endl;
//
//    // Mark neighbours of faces to be collapsed, i.e. get the first layer
//    // of triangles outside collapsedFaces.
//    // neighbours actually contains the
//    // edge with which triangle connects to collapsedFaces.
//
707
//    Map<label> neighbours;
708 709 710
//
//    labelList collapsed = collapsedFaces.toc();
//
711
//    for (const label facei : collapsed)
712
//    {
713
//        const labelList& myEdges = surf.faceEdges()[facei];
714
//
715
//        Pout<< "collapsing facei:" << facei << " uses edges:" << myEdges
716 717 718 719 720 721 722 723 724 725 726 727
//            << endl;
//
//        forAll(myEdges, myEdgeI)
//        {
//            const labelList& myFaces = surf.edgeFaces()[myEdges[myEdgeI]];
//
//            Pout<< "Edge " << myEdges[myEdgeI] << " is used by faces "
//                << myFaces << endl;
//
//            if ((myEdges[myEdgeI] != edgeI) && (myFaces.size() == 2))
//            {
//                // Get the other face
728
//                label neighbourFacei = myFaces[0];
729
//
730
//                if (neighbourFacei == facei)
731
//                {
732
//                    neighbourFacei = myFaces[1];
733 734 735
//                }
//
//                // Is 'outside' face if not in collapsedFaces itself
736
//                if (!collapsedFaces.found(neighbourFacei))
737
//                {
738
//                    neighbours.insert(neighbourFacei, myEdges[myEdgeI]);
739 740 741 742 743 744 745
//                }
//            }
//        }
//    }
//
//    // Now neighbours contains first layer of triangles outside of
//    // collapseFaces
746
//    // There should be
747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766
//    // -two if edgeI is a boundary edge
//    // since the outside 'edge' of collapseFaces should
//    // form a triangle and the face connected to edgeI is not inserted.
//    // -four if edgeI is not a boundary edge since then the outside edge forms
//    // a diamond.
//
//    // Check if any of the faces in neighbours share an edge. (n^2)
//
//    labelList neighbourList = neighbours.toc();
//
//Pout<< "edgeI:" << edgeI << "  neighbourList:" << neighbourList << endl;
//
//
//    forAll(neighbourList, i)
//    {
//        const labelList& faceIEdges = surf.faceEdges()[neighbourList[i]];
//
//        for (label j = i+1; j < neighbourList.size(); i++)
//        {
//            const labelList& faceJEdges = surf.faceEdges()[neighbourList[j]];
767
//
768
//            // Check if facei and faceJ share an edge
769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797
//            forAll(faceIEdges, fI)
//            {
//                forAll(faceJEdges, fJ)
//                {
//                    Pout<< " comparing " << faceIEdges[fI] << " to "
//                        << faceJEdges[fJ] << endl;
//
//                    if (faceIEdges[fI] == faceJEdges[fJ])
//                    {
//                        return true;
//                    }
//                }
//            }
//        }
//    }
//    Pout<< "Found no match. Returning false" << endl;
//    return false;
//}


// Finds the triangle edge cut by the plane between a point inside / on edge
// of a triangle and a point outside. Returns:
//  - cut through edge/point
//  - miss()
Foam::surfaceLocation Foam::triSurfaceTools::cutEdge
(
    const triSurface& s,
    const label triI,
    const label excludeEdgeI,
798
    const label excludePointi,
799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814

    const point& triPoint,
    const plane& cutPlane,
    const point& toPoint
)
{
    const pointField& points = s.points();
    const labelledTri& f = s[triI];
    const labelList& fEdges = s.faceEdges()[triI];

    // Get normal distance to planeN
    FixedList<scalar, 3> d;

    scalar norm = 0;
    forAll(d, fp)
    {
815
        d[fp] = cutPlane.signedDistance(points[f[fp]]);
816 817 818 819 820 821 822 823
        norm += mag(d[fp]);
    }

    // Normalise and truncate
    forAll(d, i)
    {
        d[i] /= norm;

824
        if (mag(d[i]) < 1e-6)
825 826 827 828 829 830 831 832
        {
            d[i] = 0.0;
        }
    }

    // Return information
    surfaceLocation cut;

833
    if (excludePointi != -1)
834 835 836
    {
        // Excluded point. Test only opposite edge.

837
        label fp0 = s.localFaces()[triI].find(excludePointi);
838 839 840

        if (fp0 == -1)
        {
841
            FatalErrorInFunction
842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897
                << " localF:" << s.localFaces()[triI] << abort(FatalError);
        }

        label fp1 = f.fcIndex(fp0);
        label fp2 = f.fcIndex(fp1);


        if (d[fp1] == 0.0)
        {
            cut.setHit();
            cut.setPoint(points[f[fp1]]);
            cut.elementType() = triPointRef::POINT;
            cut.setIndex(s.localFaces()[triI][fp1]);
        }
        else if (d[fp2] == 0.0)
        {
            cut.setHit();
            cut.setPoint(points[f[fp2]]);
            cut.elementType() = triPointRef::POINT;
            cut.setIndex(s.localFaces()[triI][fp2]);
        }
        else if
        (
            (d[fp1] < 0 && d[fp2] < 0)
         || (d[fp1] > 0 && d[fp2] > 0)
        )
        {
            // Both same sign. Not crossing edge at all.
            // cut already set to miss().
        }
        else
        {
            cut.setHit();
            cut.setPoint
            (
                (d[fp2]*points[f[fp1]] - d[fp1]*points[f[fp2]])
              / (d[fp2] - d[fp1])
            );
            cut.elementType() = triPointRef::EDGE;
            cut.setIndex(fEdges[fp1]);
        }
    }
    else
    {
        // Find the two intersections
        FixedList<surfaceLocation, 2> inters;
        label interI = 0;

        forAll(f, fp0)
        {
            label fp1 = f.fcIndex(fp0);

            if (d[fp0] == 0)
            {
                if (interI >= 2)
                {
898
                    FatalErrorInFunction
899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916
                        << "problem : triangle has three intersections." << nl
                        << "triangle:" << f.tri(points)
                        << " d:" << d << abort(FatalError);
                }
                inters[interI].setHit();
                inters[interI].setPoint(points[f[fp0]]);
                inters[interI].elementType() = triPointRef::POINT;
                inters[interI].setIndex(s.localFaces()[triI][fp0]);
                interI++;
            }
            else if
            (
                (d[fp0] < 0 && d[fp1] > 0)
             || (d[fp0] > 0 && d[fp1] < 0)
            )
            {
                if (interI >= 2)
                {
917
                    FatalErrorInFunction
918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017
                        << "problem : triangle has three intersections." << nl
                        << "triangle:" << f.tri(points)
                        << " d:" << d << abort(FatalError);
                }
                inters[interI].setHit();
                inters[interI].setPoint
                (
                    (d[fp0]*points[f[fp1]] - d[fp1]*points[f[fp0]])
                  / (d[fp0] - d[fp1])
                );
                inters[interI].elementType() = triPointRef::EDGE;
                inters[interI].setIndex(fEdges[fp0]);
                interI++;
            }
        }


        if (interI == 0)
        {
            // Return miss
        }
        else if (interI == 1)
        {
            // Only one intersection. Should not happen!
            cut = inters[0];
        }
        else if (interI == 2)
        {
            // Handle excludeEdgeI
            if
            (
                inters[0].elementType() == triPointRef::EDGE
             && inters[0].index() == excludeEdgeI
            )
            {
                cut = inters[1];
            }
            else if
            (
                inters[1].elementType() == triPointRef::EDGE
             && inters[1].index() == excludeEdgeI
            )
            {
                cut = inters[0];
            }
            else
            {
                // Two cuts. Find nearest.
                if
                (
                    magSqr(inters[0].rawPoint() - toPoint)
                  < magSqr(inters[1].rawPoint() - toPoint)
                )
                {
                    cut = inters[0];
                }
                else
                {
                    cut = inters[1];
                }
            }
        }
    }
    return cut;
}


// 'Snap' point on to endPoint.
void Foam::triSurfaceTools::snapToEnd
(
    const triSurface& s,
    const surfaceLocation& end,
    surfaceLocation& current
)
{
    if (end.elementType() == triPointRef::NONE)
    {
        if (current.elementType() == triPointRef::NONE)
        {
            // endpoint on triangle; current on triangle
            if (current.index() == end.index())
            {
                //if (debug)
                //{
                //    Pout<< "snapToEnd : snapping:" << current << " onto:"
                //        << end << endl;
                //}
                current = end;
                current.setHit();
            }
        }
        // No need to handle current on edge/point since tracking handles this.
    }
    else if (end.elementType() == triPointRef::EDGE)
    {
        if (current.elementType() == triPointRef::NONE)
        {
            // endpoint on edge; current on triangle
            const labelList& fEdges = s.faceEdges()[current.index()];

1018
            if (fEdges.found(end.index()))
1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064
            {
                //if (debug)
                //{
                //    Pout<< "snapToEnd : snapping:" << current << " onto:"
                //        << end << endl;
                //}
                current = end;
                current.setHit();
            }
        }
        else if (current.elementType() == triPointRef::EDGE)
        {
            // endpoint on edge; current on edge
            if (current.index() == end.index())
            {
                //if (debug)
                //{
                //    Pout<< "snapToEnd : snapping:" << current << " onto:"
                //        << end << endl;
                //}
                current = end;
                current.setHit();
            }
        }
        else
        {
            // endpoint on edge; current on point
            const edge& e = s.edges()[end.index()];

            if (current.index() == e[0] || current.index() == e[1])
            {
                //if (debug)
                //{
                //    Pout<< "snapToEnd : snapping:" << current << " onto:"
                //        << end << endl;
                //}
                current = end;
                current.setHit();
            }
        }
    }
    else    // end.elementType() == POINT
    {
        if (current.elementType() == triPointRef::NONE)
        {
            // endpoint on point; current on triangle
1065
            const triSurface::face_type& f = s.localFaces()[current.index()];
1066

1067
            if (f.found(end.index()))
1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121
            {
                //if (debug)
                //{
                //    Pout<< "snapToEnd : snapping:" << current << " onto:"
                //        << end << endl;
                //}
                current = end;
                current.setHit();
            }
        }
        else if (current.elementType() == triPointRef::EDGE)
        {
            // endpoint on point; current on edge
            const edge& e = s.edges()[current.index()];

            if (end.index() == e[0] || end.index() == e[1])
            {
                //if (debug)
                //{
                //    Pout<< "snapToEnd : snapping:" << current << " onto:"
                //        << end << endl;
                //}
                current = end;
                current.setHit();
            }
        }
        else
        {
            // endpoint on point; current on point
            if (current.index() == end.index())
            {
                //if (debug)
                //{
                //    Pout<< "snapToEnd : snapping:" << current << " onto:"
                //        << end << endl;
                //}
                current = end;
                current.setHit();
            }
        }
    }
}


// Start:
//  - location
//  - element type (triangle/edge/point) and index
//  - triangle to exclude
Foam::surfaceLocation Foam::triSurfaceTools::visitFaces
(
    const triSurface& s,
    const labelList& eFaces,
    const surfaceLocation& start,
    const label excludeEdgeI,
1122
    const label excludePointi,
1123 1124 1125 1126 1127 1128 1129 1130
    const surfaceLocation& end,
    const plane& cutPlane
)
{
    surfaceLocation nearest;

    scalar minDistSqr = Foam::sqr(GREAT);

1131
    for (const label triI : eFaces)
1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142
    {
        // Make sure we don't revisit previous face
        if (triI != start.triangle())
        {
            if (end.elementType() == triPointRef::NONE && end.index() == triI)
            {
                // Endpoint is in this triangle. Jump there.
                nearest = end;
                nearest.setHit();
                nearest.triangle() = triI;
                break;
1143
            }
1144 1145 1146 1147 1148 1149 1150 1151 1152
            else
            {
               // Which edge is cut.

                surfaceLocation cutInfo = cutEdge
                (
                    s,
                    triI,
                    excludeEdgeI,       // excludeEdgeI
1153
                    excludePointi,      // excludePointi
1154 1155 1156 1157 1158 1159 1160 1161
                    start.rawPoint(),
                    cutPlane,
                    end.rawPoint()
                );

                // If crossing an edge we expect next edge to be cut.
                if (excludeEdgeI != -1 && !cutInfo.hit())
                {
1162
                    FatalErrorInFunction
1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188
                        << "Triangle:" << triI
                        << " excludeEdge:" << excludeEdgeI
                        << " point:" << start.rawPoint()
                        << " plane:" << cutPlane
                        << " . No intersection!" << abort(FatalError);
                }

                if (cutInfo.hit())
                {
                    scalar distSqr = magSqr(cutInfo.rawPoint()-end.rawPoint());

                    if (distSqr < minDistSqr)
                    {
                        minDistSqr = distSqr;
                        nearest = cutInfo;
                        nearest.triangle() = triI;
                        nearest.setMiss();
                    }
                }
            }
        }
    }

    if (nearest.triangle() == -1)
    {
        // Did not move from edge. Give warning? Return something special?
Andrew Heather's avatar
Andrew Heather committed
1189
        // For now responsibility of caller to make sure that nothing has
1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208
        // moved.
    }

    return nearest;
}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //


// Write pointField to file
void Foam::triSurfaceTools::writeOBJ
(
    const fileName& fName,
    const pointField& pts
)
{
    OFstream outFile(fName);

1209
    for (const point& pt : pts)
1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272
    {
        outFile<< "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
    }
    Pout<< "Written " << pts.size() << " vertices to file " << fName << endl;
}


// Write vertex subset to OBJ format file
void Foam::triSurfaceTools::writeOBJ
(
    const triSurface& surf,
    const fileName& fName,
    const boolList& markedVerts
)
{
    OFstream outFile(fName);

    label nVerts = 0;
    forAll(markedVerts, vertI)
    {
        if (markedVerts[vertI])
        {
            const point& pt = surf.localPoints()[vertI];

            outFile<< "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;

            nVerts++;
        }
    }
    Pout<< "Written " << nVerts << " vertices to file " << fName << endl;
}


// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// Addressing helper functions:


// Get all triangles using vertices of edge
void Foam::triSurfaceTools::getVertexTriangles
(
    const triSurface& surf,
    const label edgeI,
    labelList& edgeTris
)
{
    const edge& e = surf.edges()[edgeI];
    const labelList& myFaces = surf.edgeFaces()[edgeI];

    label face1I = myFaces[0];
    label face2I = -1;
    if (myFaces.size() == 2)
    {
        face2I = myFaces[1];
    }

    const labelList& startFaces = surf.pointFaces()[e.start()];
    const labelList& endFaces = surf.pointFaces()[e.end()];

    // Number of triangles is sum of pointfaces - common faces
    // (= faces using edge)
    edgeTris.setSize(startFaces.size() + endFaces.size() - myFaces.size());

    label nTris = 0;
1273
    for (const label facei : startFaces)
1274
    {
1275
        edgeTris[nTris++] = facei;
1276 1277
    }

1278
    for (const label facei : endFaces)
1279
    {
1280
        if ((facei != face1I) && (facei != face2I))
1281
        {
1282
            edgeTris[nTris++] = facei;
1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295
        }
    }
}


// Get all vertices connected to vertices of edge
Foam::labelList Foam::triSurfaceTools::getVertexVertices
(
    const triSurface& surf,
    const edge& e
)
{
    const edgeList& edges = surf.edges();
1296 1297
    const label v1 = e.start();
    const label v2 = e.end();
1298 1299 1300 1301 1302 1303

    // Get all vertices connected to v1 or v2 through an edge
    labelHashSet vertexNeighbours;

    const labelList& v1Edges = surf.pointEdges()[v1];

1304
    for (const label edgei : v1Edges)
1305
    {
1306
        vertexNeighbours.insert(edges[edgei].otherVertex(v1));
1307 1308 1309 1310
    }

    const labelList& v2Edges = surf.pointEdges()[v2];

1311
    for (const label edgei : v2Edges)
1312
    {
1313
        vertexNeighbours.insert(edges[edgei].otherVertex(v2));
1314 1315 1316 1317 1318 1319 1320 1321 1322
    }
    return vertexNeighbours.toc();
}


// Get the other face using edgeI
Foam::label Foam::triSurfaceTools::otherFace
(
    const triSurface& surf,
1323
    const label facei,
1324 1325 1326 1327 1328 1329 1330 1331 1332
    const label edgeI
)
{
    const labelList& myFaces = surf.edgeFaces()[edgeI];

    if (myFaces.size() != 2)
    {
        return -1;
    }
1333 1334 1335 1336
    else if (facei == myFaces[0])
    {
        return myFaces[1];
    }
1337 1338
    else
    {
1339
        return myFaces[0];
1340 1341 1342 1343
    }
}


1344
// Get the two edges on facei counterclockwise after edgeI
1345 1346 1347
void Foam::triSurfaceTools::otherEdges
(
    const triSurface& surf,
1348
    const label facei,
1349 1350 1351 1352 1353
    const label edgeI,
    label& e1,
    label& e2
)
{
1354
    const labelList& eFaces = surf.faceEdges()[facei];
1355

1356
    label i0 = eFaces.find(edgeI);
1357 1358 1359

    if (i0 == -1)
    {
1360 1361
        FatalErrorInFunction
            << "Edge " << surf.edges()[edgeI] << " not in face "
1362
            << surf.localFaces()[facei] << abort(FatalError);
1363 1364 1365 1366 1367 1368 1369 1370 1371 1372
    }

    label i1 = eFaces.fcIndex(i0);
    label i2 = eFaces.fcIndex(i1);

    e1 = eFaces[i1];
    e2 = eFaces[i2];
}


1373
// Get the two vertices on facei counterclockwise vertI
1374 1375 1376
void Foam::triSurfaceTools::otherVertices
(
    const triSurface& surf,
1377
    const label facei,
1378 1379 1380 1381 1382
    const label vertI,
    label& vert1I,
    label& vert2I
)
{
1383
    const labelledTri& f = surf.localFaces()[facei];
1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401

    if (vertI == f[0])
    {
        vert1I = f[1];
        vert2I = f[2];
    }
    else if (vertI == f[1])
    {
        vert1I = f[2];
        vert2I = f[0];
    }
    else if (vertI == f[2])
    {
        vert1I = f[0];
        vert2I = f[1];
    }
    else
    {
1402
        FatalErrorInFunction
1403 1404
            << "Vertex " << vertI << " not in face " << f << nl
            << abort(FatalError);
1405 1406 1407 1408 1409 1410 1411 1412
    }
}


// Get edge opposite vertex
Foam::label Foam::triSurfaceTools::oppositeEdge
(
    const triSurface& surf,
1413
    const label facei,
1414 1415 1416
    const label vertI
)
{
1417
    const labelList& myEdges = surf.faceEdges()[facei];
1418

1419
    for (const label edgei : myEdges)
1420
    {
1421
        const edge& e = surf.edges()[edgei];
1422