isoSurfaceTopo.C 40.5 KB
Newer Older
1 2 3
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4
   \\    /   O peration     |
OpenFOAM bot's avatar
OpenFOAM bot committed
5
    \\  /    A nd           | www.openfoam.com
6
     \\/     M anipulation  |
OpenFOAM bot's avatar
OpenFOAM bot committed
7
-------------------------------------------------------------------------------
OpenFOAM bot's avatar
OpenFOAM bot committed
8 9
    Copyright (C) 2011-2019 OpenFOAM Foundation
    Copyright (C) 2019 OpenCFD Ltd.
10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    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.

    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
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

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

#include "isoSurfaceTopo.H"
30
#include "isoSurface.H"
31 32 33 34
#include "polyMesh.H"
#include "tetMatcher.H"
#include "tetPointRef.H"
#include "DynamicField.H"
35 36
#include "syncTools.H"
#include "polyMeshTetDecomposition.H"
37
#include "cyclicACMIPolyPatch.H"
38 39 40 41 42 43 44 45 46

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

namespace Foam
{
    defineTypeNameAndDebug(isoSurfaceTopo, 0);
}


47
// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
48

49
namespace Foam
50
{
51 52 53 54 55 56 57 58 59 60 61 62 63
    // Does any edge of triangle cross iso value?
    inline static bool isTriCut
    (
        const label a,
        const label b,
        const label c,
        const scalarField& pointValues,
        const scalar isoValue
    )
    {
        const bool aLower = (pointValues[a] < isoValue);
        const bool bLower = (pointValues[b] < isoValue);
        const bool cLower = (pointValues[c] < isoValue);
64

65 66
        return !(aLower == bLower && aLower == cLower);
    }
67 68 69
}


70 71
// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //

72 73 74 75 76 77
Foam::isoSurfaceTopo::cellCutType Foam::isoSurfaceTopo::calcCutType
(
    const bool isTet,
    const label celli
) const
{
78 79 80 81 82
    if (ignoreCells_.test(celli))
    {
        return NOTCUT;
    }

83 84 85 86 87 88
    const cell& cFaces = mesh_.cells()[celli];

    if (isTet)
    {
        for (const label facei : cFaces)
        {
89 90 91 92 93 94 95 96 97
            if
            (
               !mesh_.isInternalFace(facei)
             && ignoreBoundaryFaces_.test(facei-mesh_.nInternalFaces())
            )
            {
                continue;
            }

98 99 100 101
            const face& f = mesh_.faces()[facei];

            for (label fp = 1; fp < f.size() - 1; ++fp)
            {
102
                if (isTriCut(f[0], f[fp], f[f.fcIndex(fp)], pVals_, iso_))
103 104 105 106 107 108 109 110 111
                {
                    return CUT;
                }
            }
        }
        return NOTCUT;
    }


112
    const bool cellLower = (cVals_[celli] < iso_);
113

114 115
    // First check if there is any cut in cell
    bool edgeCut = false;
116

117 118
    for (const label facei : cFaces)
    {
119 120 121 122 123 124 125 126 127
        if
        (
           !mesh_.isInternalFace(facei)
         && ignoreBoundaryFaces_.test(facei-mesh_.nInternalFaces())
        )
        {
            continue;
        }

128
        const face& f = mesh_.faces()[facei];
129

130 131 132 133
        // Check pyramids cut
        for (const label pointi : f)
        {
            if ((pVals_[pointi] < iso_) != cellLower)
134
            {
135 136
                edgeCut = true;
                break;
137
            }
138
        }
139

140 141 142 143
        if (edgeCut)
        {
            break;
        }
144

145 146 147
        // Check (triangulated) face edges
        // With fallback for problem decompositions
        const label fp0 = (tetBasePtIs_[facei] < 0 ? 0 : tetBasePtIs_[facei]);
148

149 150 151 152
        label fp = f.fcIndex(fp0);
        for (label i = 2; i < f.size(); ++i)
        {
            const label nextFp = f.fcIndex(fp);
153

154
            if (isTriCut(f[fp0], f[fp], f[nextFp], pVals_, iso_))
155
            {
156
                edgeCut = true;
157 158
                break;
            }
159 160

            fp = nextFp;
161 162 163 164
        }

        if (edgeCut)
        {
165 166 167
            break;
        }
    }
168

169 170 171 172 173
    if (edgeCut)
    {
        // Count actual cuts (expensive since addressing needed)
        // Note: not needed if you don't want to preserve maxima/minima
        // centred around cellcentre. In that case just always return CUT
174

175
        const labelList& cPoints = mesh_.cellPoints(celli);
176

177
        label nPyrCuts = 0;
178

179 180 181
        for (const label pointi : cPoints)
        {
            if ((pVals_[pointi] < iso_) != cellLower)
182
            {
183
                ++nPyrCuts;
184 185
            }
        }
186 187 188 189 190 191 192 193 194

        if (nPyrCuts == cPoints.size())
        {
            return SPHERE;
        }
        else if (nPyrCuts)
        {
            return CUT;
        }
195
    }
196 197

    return NOTCUT;
198 199 200 201 202 203 204 205 206
}


Foam::label Foam::isoSurfaceTopo::calcCutTypes
(
    tetMatcher& tet,
    List<cellCutType>& cellCutTypes
)
{
207
    cellCutTypes.setSize(mesh_.nCells());
208 209

    label nCutCells = 0;
210
    forAll(cellCutTypes, celli)
211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228
    {
        cellCutTypes[celli] = calcCutType(tet.isA(mesh_, celli), celli);

        if (cellCutTypes[celli] == CUT)
        {
            ++nCutCells;
        }
    }

    if (debug)
    {
        Pout<< "isoSurfaceCell : candidate cut cells "
            << nCutCells << " / " << mesh_.nCells() << endl;
    }
    return nCutCells;
}


229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276
Foam::scalar Foam::isoSurfaceTopo::minTetQ
(
    const label facei,
    const label faceBasePtI
) const
{
    scalar q = polyMeshTetDecomposition::minQuality
    (
        mesh_,
        mesh_.cellCentres()[mesh_.faceOwner()[facei]],
        facei,
        true,
        faceBasePtI
    );

    if (mesh_.isInternalFace(facei))
    {
        q = min
        (
            q,
            polyMeshTetDecomposition::minQuality
            (
                mesh_,
                mesh_.cellCentres()[mesh_.faceNeighbour()[facei]],
                facei,
                false,
                faceBasePtI
            )
        );
    }
    return q;
}


void Foam::isoSurfaceTopo::fixTetBasePtIs()
{
    // Determine points used by two faces on the same cell
    const cellList& cells = mesh_.cells();
    const faceList& faces = mesh_.faces();
    const labelList& faceOwner = mesh_.faceOwner();
    const labelList& faceNeighbour = mesh_.faceNeighbour();


    // Get face triangulation base point
    tetBasePtIs_ = mesh_.tetBasePtIs();


    // Pre-filter: mark all cells with illegal base points
277 278
    bitSet problemCells(cells.size());

279 280 281 282
    forAll(tetBasePtIs_, facei)
    {
        if (tetBasePtIs_[facei] == -1)
        {
283 284
            problemCells.set(faceOwner[facei]);

285 286
            if (mesh_.isInternalFace(facei))
            {
287
                problemCells.set(faceNeighbour[facei]);
288 289 290 291 292
            }
        }
    }


293 294 295
    // Mark all points that are shared by just two faces within an adjacent
    // problem cell as problematic
    bitSet problemPoints(mesh_.points().size());
296 297

    {
298 299 300
        // Number of times a point occurs in a cell.
        // Used to detect dangling vertices (count = 2)
        Map<label> pointCount;
301

302 303 304
        // Analyse problem cells for points shared by two faces only
        for (const label celli : problemCells)
        {
305 306
            pointCount.clear();

307
            for (const label facei : cells[celli])
308
            {
309
                for (const label pointi : faces[facei])
310
                {
311
                    ++pointCount(pointi);
312 313 314
                }
            }

315
            forAllConstIters(pointCount, iter)
316
            {
317
                if (iter.val() == 1)
318
                {
319 320
                    FatalErrorInFunction
                        << "point:" << iter.key()
321
                        << " at:" << mesh_.points()[iter.key()]
322 323
                        << " only used by one face" << nl
                        << exit(FatalError);
324
                }
325
                else if (iter.val() == 2)
326
                {
327
                    problemPoints.set(iter.key());
328 329
                }
            }
330 331
        }
    }
332

333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365

    // For all faces which form a part of a problem-cell, check if the base
    // point is adjacent to any problem points. If it is, re-calculate the base
    // point so that it is not.
    label nAdapted = 0;
    forAll(tetBasePtIs_, facei)
    {
        if
        (
            problemCells[faceOwner[facei]]
         || (mesh_.isInternalFace(facei) && problemCells[faceNeighbour[facei]])
        )
        {
            const face& f = faces[facei];

            // Check if either of the points adjacent to the base point is a
            // problem point. If not, the existing base point can be retained.
            const label fp0 = tetBasePtIs_[facei] < 0 ? 0 : tetBasePtIs_[facei];

            if
            (
                !problemPoints[f[f.rcIndex(fp0)]]
             && !problemPoints[f[f.fcIndex(fp0)]]
            )
            {
                continue;
            }

            // A new base point is required. Pick the point that results in the
            // least-worst tet and which is not adjacent to any problem points.
            scalar maxQ = -GREAT;
            label maxFp = -1;
            forAll(f, fp)
366
            {
367 368 369 370 371
                if
                (
                    !problemPoints[f[f.rcIndex(fp)]]
                 && !problemPoints[f[f.fcIndex(fp)]]
                )
372
                {
373 374
                    const scalar q = minTetQ(facei, fp);
                    if (q > maxQ)
375
                    {
376 377
                        maxQ = q;
                        maxFp = fp;
378 379 380
                    }
                }
            }
381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396

            if (maxFp != -1)
            {
                // Success! Set the new base point
                tetBasePtIs_[facei] = maxFp;
            }
            else
            {
                // No point was found on face that would not result in some
                // duplicate triangle. Do what? Continue and hope? Spit an
                // error? Silently or noisily reduce the filtering level?

                tetBasePtIs_[facei] = 0;
            }

            ++nAdapted;
397 398 399
        }
    }

400

401 402 403 404 405 406 407 408 409 410
    if (debug)
    {
        Pout<< "isoSurface : adapted starting point of triangulation on "
            << nAdapted << " faces." << endl;
    }

    syncTools::syncFaceList(mesh_, tetBasePtIs_, maxEqOp<label>());
}


411 412 413 414 415 416 417 418 419 420 421 422
Foam::label Foam::isoSurfaceTopo::generatePoint
(
    const label facei,
    const bool edgeIsDiag,
    const edge& vertices,

    DynamicList<edge>& pointToVerts,
    DynamicList<label>& pointToFace,
    DynamicList<bool>& pointFromDiag,
    EdgeMap<label>& vertsToPoint
) const
{
423 424
    const auto edgeFnd = vertsToPoint.cfind(vertices);
    if (edgeFnd.found())
425
    {
426
        return edgeFnd.val();
427
    }
428 429 430 431 432 433 434 435 436

    // Generate new point
    label pointi = pointToVerts.size();

    pointToVerts.append(vertices);
    pointToFace.append(facei);
    pointFromDiag.append(edgeIsDiag);
    vertsToPoint.insert(vertices, pointi);
    return pointi;
437 438 439 440 441 442 443 444 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 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 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 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 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 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 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 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845
}


void Foam::isoSurfaceTopo::generateTriPoints
(
    const label facei,
    const FixedList<scalar, 4>& s,
    const FixedList<point, 4>& p,
    const FixedList<label, 4>& pIndex,
    const FixedList<bool, 6>& edgeIsDiag,// Per tet edge whether is face diag

    DynamicList<edge>& pointToVerts,
    DynamicList<label>& pointToFace,
    DynamicList<bool>& pointFromDiag,

    EdgeMap<label>& vertsToPoint,
    DynamicList<label>& verts       // Every three verts is new triangle
) const
{
    int triIndex = 0;
    if (s[0] < iso_)
    {
        triIndex |= 1;
    }
    if (s[1] < iso_)
    {
        triIndex |= 2;
    }
    if (s[2] < iso_)
    {
        triIndex |= 4;
    }
    if (s[3] < iso_)
    {
        triIndex |= 8;
    }

    // Form the vertices of the triangles for each case
    switch (triIndex)
    {
        case 0x00:
        case 0x0F:
        break;

        case 0x01:
        case 0x0E:
        {
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[0],
                    edge(pIndex[0], pIndex[1]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[1],
                    edge(pIndex[0], pIndex[2]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[2],
                    edge(pIndex[0], pIndex[3]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );

            if (triIndex == 0x0E)
            {
                // Flip normals
                label sz = verts.size();
                Swap(verts[sz-2], verts[sz-1]);
            }
        }
        break;

        case 0x02:
        case 0x0D:
        {
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[0],
                    edge(pIndex[1], pIndex[0]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[3],
                    edge(pIndex[1], pIndex[3]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[4],
                    edge(pIndex[1], pIndex[2]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );

            if (triIndex == 0x0D)
            {
                // Flip normals
                label sz = verts.size();
                Swap(verts[sz-2], verts[sz-1]);
            }
        }
        break;

        case 0x03:
        case 0x0C:
        {
            label p0p2
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[1],
                    edge(pIndex[0], pIndex[2]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            label p1p3
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[3],
                    edge(pIndex[1], pIndex[3]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );

            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[2],
                    edge(pIndex[0], pIndex[3]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append(p1p3);
            verts.append(p0p2);
            verts.append(p1p3);
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[4],
                    edge(pIndex[1], pIndex[2]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append(p0p2);

            if (triIndex == 0x0C)
            {
                // Flip normals
                label sz = verts.size();
                Swap(verts[sz-5], verts[sz-4]);
                Swap(verts[sz-2], verts[sz-1]);
            }
        }
        break;

        case 0x04:
        case 0x0B:
        {
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[1],
                    edge(pIndex[2], pIndex[0]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[4],
                    edge(pIndex[2], pIndex[1]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[5],
                    edge(pIndex[2], pIndex[3]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );

            if (triIndex == 0x0B)
            {
                // Flip normals
                label sz = verts.size();
                Swap(verts[sz-2], verts[sz-1]);
            }
        }
        break;

        case 0x05:
        case 0x0A:
        {
            label p0p1
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[0],
                    edge(pIndex[0], pIndex[1]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            label p2p3
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[5],
                    edge(pIndex[2], pIndex[3]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );

            verts.append(p0p1);
            verts.append(p2p3);
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[2],
                    edge(pIndex[0], pIndex[3]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append(p0p1);
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[4],
                    edge(pIndex[1], pIndex[2]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append(p2p3);

            if (triIndex == 0x0A)
            {
                // Flip normals
                label sz = verts.size();
                Swap(verts[sz-5], verts[sz-4]);
                Swap(verts[sz-2], verts[sz-1]);
            }
        }
        break;

        case 0x06:
        case 0x09:
        {
            label p0p1
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[0],
                    edge(pIndex[0], pIndex[1]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            label p2p3
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[5],
                    edge(pIndex[2], pIndex[3]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );

            verts.append(p0p1);
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[3],
                    edge(pIndex[1], pIndex[3]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append(p2p3);
            verts.append(p0p1);
            verts.append(p2p3);
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[1],
                    edge(pIndex[0], pIndex[2]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );

            if (triIndex == 0x09)
            {
                // Flip normals
                label sz = verts.size();
                Swap(verts[sz-5], verts[sz-4]);
                Swap(verts[sz-2], verts[sz-1]);
            }
        }
        break;

        case 0x08:
        case 0x07:
        {
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[2],
                    edge(pIndex[3], pIndex[0]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[5],
                    edge(pIndex[3], pIndex[2]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            verts.append
            (
                generatePoint
                (
                    facei,
                    edgeIsDiag[3],
                    edge(pIndex[3], pIndex[1]),
                    pointToVerts, pointToFace, pointFromDiag, vertsToPoint
                )
            );
            if (triIndex == 0x07)
            {
                // Flip normals
                label sz = verts.size();
                Swap(verts[sz-2], verts[sz-1]);
            }
        }
        break;
    }
}


void Foam::isoSurfaceTopo::generateTriPoints
(
    const label celli,
    const bool isTet,

    DynamicList<edge>& pointToVerts,
    DynamicList<label>& pointToFace,
    DynamicList<bool>& pointFromDiag,

    EdgeMap<label>& vertsToPoint,
    DynamicList<label>& verts,
    DynamicList<label>& faceLabels
) const
{
846 847 848 849
    const faceList& faces = mesh_.faces();
    const labelList& faceOwner = mesh_.faceOwner();
    const pointField& points = mesh_.points();
    const cell& cFaces = mesh_.cells()[celli];
850 851 852 853 854 855 856

    if (isTet)
    {
        // For tets don't do cell-centre decomposition, just use the
        // tet points and values

        label facei = cFaces[0];
857
        const face& f0 = faces[facei];
858

859 860
        // Get the other point from f1. Tbd: check if not duplicate face
        // (ACMI / ignoreBoundaryFaces_).
861
        const face& f1 = faces[cFaces[1]];
862 863 864 865
        label oppositeI = -1;
        forAll(f1, fp)
        {
            oppositeI = f1[fp];
866
            if (!f0.found(oppositeI))
867 868 869 870 871 872 873 874 875 876 877
            {
                break;
            }
        }


        label p0 = f0[0];
        label p1 = f0[1];
        label p2 = f0[2];
        FixedList<bool, 6> edgeIsDiag(false);

878
        if (faceOwner[facei] == celli)
879 880 881 882 883 884
        {
            Swap(p1, p2);
        }

        tetPointRef tet
        (
885 886 887 888
            points[p0],
            points[p1],
            points[p2],
            points[oppositeI]
889 890 891 892 893 894 895 896 897 898 899 900 901 902 903
        );

        label startTrii = verts.size();
        generateTriPoints
        (
            facei,
            FixedList<scalar, 4>
            ({
                pVals_[p0],
                pVals_[p1],
                pVals_[p2],
                pVals_[oppositeI]
            }),
            FixedList<point, 4>
            ({
904 905 906 907
                points[p0],
                points[p1],
                points[p2],
                points[oppositeI]
908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932
            }),
            FixedList<label, 4>
            ({
                p0,
                p1,
                p2,
                oppositeI
            }),
            edgeIsDiag,

            pointToVerts,
            pointToFace,
            pointFromDiag,
            vertsToPoint,
            verts       // Every three verts is new triangle
        );

        label nTris = (verts.size()-startTrii)/3;
        for (label i = 0; i < nTris; ++i)
        {
            faceLabels.append(facei);
        }
    }
    else
    {
933 934
        const pointField& cellCentres = mesh_.cellCentres();

935 936
        for (const label facei : cFaces)
        {
937 938 939 940 941 942 943 944 945
            if
            (
               !mesh_.isInternalFace(facei)
             && ignoreBoundaryFaces_.test(facei-mesh_.nInternalFaces())
            )
            {
                continue;
            }

946
            const face& f = faces[facei];
947

948
            label fp0 = tetBasePtIs_[facei];
949 950 951

            label startTrii = verts.size();

952
            // Fallback
953 954 955 956 957 958 959 960 961 962 963 964 965 966 967
            if (fp0 < 0)
            {
                fp0 = 0;
            }

            label fp = f.fcIndex(fp0);
            for (label i = 2; i < f.size(); ++i)
            {
                label nextFp = f.fcIndex(fp);

                FixedList<bool, 6> edgeIsDiag(false);

                label p0 = f[fp0];
                label p1 = f[fp];
                label p2 = f[nextFp];
968
                if (faceOwner[facei] == celli)
969 970 971 972 973 974 975 976 977 978 979 980 981
                {
                    Swap(p1, p2);
                    if (i != 2) edgeIsDiag[1] = true;
                    if (i != f.size()-1) edgeIsDiag[0] = true;
                }
                else
                {
                    if (i != 2) edgeIsDiag[0] = true;
                    if (i != f.size()-1) edgeIsDiag[1] = true;
                }

                tetPointRef tet
                (
982 983 984 985
                    points[p0],
                    points[p1],
                    points[p2],
                    cellCentres[celli]
986 987 988 989 990 991 992 993 994 995 996 997 998 999
                );

                generateTriPoints
                (
                    facei,
                    FixedList<scalar, 4>
                    ({
                        pVals_[p0],
                        pVals_[p1],
                        pVals_[p2],
                        cVals_[celli]
                    }),
                    FixedList<point, 4>
                    ({
1000 1001 1002 1003
                        points[p0],
                        points[p1],
                        points[p2],
                        cellCentres[celli]
1004 1005 1006 1007 1008 1009
                    }),
                    FixedList<label, 4>
                    ({
                        p0,
                        p1,
                        p2,
1010
                        mesh_.nPoints()+celli
1011 1012 1013 1014 1015 1016 1017 1018 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 1065 1066 1067 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 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 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 1189 1190 1191 1192 1193 1194 1195 1196
                    }),
                    edgeIsDiag,

                    pointToVerts,
                    pointToFace,
                    pointFromDiag,
                    vertsToPoint,
                    verts       // Every three verts is new triangle
                );

                fp = nextFp;
            }

            label nTris = (verts.size()-startTrii)/3;
            for (label i = 0; i < nTris; ++i)
            {
                faceLabels.append(facei);
            }
        }
    }
}


void Foam::isoSurfaceTopo::triangulateOutside
(
    const bool filterDiag,
    const primitivePatch& pp,
    const boolList& pointFromDiag,
    const labelList& pointToFace,
    const label cellID,

    DynamicList<face>& compactFaces,
    DynamicList<label>& compactCellIDs
) const
{
    // We can form pockets:
    // - 1. triangle on face
    // - 2. multiple triangles on interior (from diag edges)
    // - the edge loop will be pocket since it is only the diag
    //   edges that give it volume?

    // Retriangulate the exterior loops
    const labelListList& edgeLoops = pp.edgeLoops();
    const labelList& mp = pp.meshPoints();

    for (const labelList& loop : edgeLoops)
    {
        if (loop.size() > 2)
        {
            compactFaces.append(face(0));
            face& f = compactFaces.last();

            f.setSize(loop.size());
            label fpi = 0;
            forAll(f, i)
            {
                label pointi = mp[loop[i]];
                if (filterDiag && pointFromDiag[pointi])
                {
                    label prevPointi = mp[loop[loop.fcIndex(i)]];
                    if
                    (
                        pointFromDiag[prevPointi]
                     && (pointToFace[pointi] != pointToFace[prevPointi])
                    )
                    {
                        f[fpi++] = pointi;
                    }
                    else
                    {
                        // Filter out diagonal point
                    }
                }
                else
                {
                    f[fpi++] = pointi;
                }
            }

            if (fpi > 2)
            {
                f.setSize(fpi);
            }
            else
            {
                // Keep original face
                forAll(f, i)
                {
                    label pointi = mp[loop[i]];
                    f[i] = pointi;
                }
            }
            compactCellIDs.append(cellID);
        }
    }
}


Foam::MeshedSurface<Foam::face> Foam::isoSurfaceTopo::removeInsidePoints
(
    const bool filterDiag,
    const MeshStorage& s,
    const boolList& pointFromDiag,
    const labelList& pointToFace,
    const labelList& start,                 // Per cell the starting triangle
    DynamicList<label>& pointCompactMap,    // Per returned point the original
    DynamicList<label>& compactCellIDs      // Per returned tri the cellID
) const
{
    const pointField& points = s.points();

    pointCompactMap.clear();
    compactCellIDs.clear();

    DynamicList<face> compactFaces(s.size()/8);

    for (label celli = 0; celli < start.size()-1; ++celli)
    {
        // All triangles of the current cell

        label nTris = start[celli+1]-start[celli];

        if (nTris)
        {
            const SubList<face> cellFaces(s, nTris, start[celli]);
            const primitivePatch pp(cellFaces, points);

            triangulateOutside
            (
                filterDiag,
                pp,
                pointFromDiag,
                pointToFace,
                //protectedFace,
                celli,

                compactFaces,
                compactCellIDs
            );
        }
    }


    // Compact out unused points
    // Pick up the used vertices
    labelList oldToCompact(points.size(), -1);
    DynamicField<point> compactPoints(points.size());
    pointCompactMap.clear();

    for (face& f : compactFaces)
    {
        forAll(f, fp)
        {
            label pointi = f[fp];
            label compacti = oldToCompact[pointi];
            if (compacti == -1)
            {
                compacti = compactPoints.size();
                oldToCompact[pointi] = compacti;
                compactPoints.append(points[pointi]);
                pointCompactMap.append(pointi);
            }
            f[fp] = compacti;
        }
    }


    MeshStorage cpSurface
    (
        std::move(compactPoints),
        std::move(compactFaces),
        s.surfZones()
    );

    return cpSurface;
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

Foam::isoSurfaceTopo::isoSurfaceTopo
(
    const polyMesh& mesh,
    const scalarField& cVals,
    const scalarField& pVals,
    const scalar iso,
1197 1198 1199
    const filterType filter,
    const boundBox& bounds,
    const bitSet& ignoreCells
1200 1201
)
:
1202
    isoSurfaceBase(iso, bounds),
1203 1204 1205
    mesh_(mesh),
    cVals_(cVals),
    pVals_(pVals),
1206
    ignoreCells_(ignoreCells)
1207 1208 1209
{
    if (debug)
    {
1210 1211 1212 1213 1214 1215 1216 1217
        Pout<< "isoSurfaceTopo::"
            << "    cell min/max  : "
            << min(cVals_) << " / "
            << max(cVals_) << nl
            << "    point min/max : "
            << min(pVals_) << " / "
            << max(pVals_) << nl
            << "    isoValue      : " << iso << nl
1218 1219
            << "    filter        : " << isoSurfaceTopo::filterNames[filter]
            << nl
1220 1221 1222 1223
            << "    mesh span     : " << mesh.bounds().mag() << nl
            << "    ignoreCells   : " << ignoreCells_.count()
            << " / " << cVals_.size() << nl
            << endl;
1224 1225
    }

1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246
    // Determine boundary pyramids to ignore (since originating from ACMI faces)
    // Maybe needs to be optional argument or more general detect colocated
    // faces.
    {
        const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
        forAll(pbm, patchi)
        {
            const polyPatch& pp = pbm[patchi];
            if (isA<cyclicACMIPolyPatch>(pp))
            {
                ignoreBoundaryFaces_.setSize(mesh_.nBoundaryFaces());
                forAll(pp, i)
                {
                    label facei = pp.start()+i;
                    ignoreBoundaryFaces_.set(facei-mesh_.nInternalFaces());
                }
            }
        }
    }


1247 1248
    fixTetBasePtIs();

1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272
    tetMatcher tet;

    // Determine if any cut through cell
    List<cellCutType> cellCutTypes;
    const label nCutCells = calcCutTypes(tet, cellCutTypes);

    // Per cell: 5 pyramids cut, each generating 2 triangles
    //  - pointToVerts : from generated iso point to originating mesh verts
    DynamicList<edge> pointToVerts(10*nCutCells);
    //  - pointToFace : from generated iso point to originating mesh face
    DynamicList<label> pointToFace(10*nCutCells);
    //  - pointToFace : from generated iso point whether is on face diagonal
    DynamicList<bool> pointFromDiag(10*nCutCells);

    // Per cell: number of intersected edges:
    //          - four faces cut so 4 mesh edges + 4 face-diagonal edges
    //          - 4 of the pyramid edges
    EdgeMap<label> vertsToPoint(12*nCutCells);
    DynamicList<label> verts(12*nCutCells);
    // Per cell: 5 pyramids cut (since only one pyramid not cut)
    DynamicList<label> faceLabels(5*nCutCells);
    DynamicList<label> cellLabels(5*nCutCells);


1273
    labelList startTri(mesh_.nCells()+1, Zero);
1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328

    for (label celli = 0; celli < mesh_.nCells(); ++celli)
    {
        startTri[celli] = faceLabels.size();
        if (cellCutTypes[celli] != NOTCUT)
        {
            generateTriPoints
            (
                celli,
                tet.isA(mesh_, celli),

                pointToVerts,
                pointToFace,
                pointFromDiag,

                vertsToPoint,
                verts,
                faceLabels
            );

            for (label i = startTri[celli]; i < faceLabels.size(); ++i)
            {
                cellLabels.append(celli);
            }
        }
    }
    startTri[mesh_.nCells()] = faceLabels.size();


    pointToVerts_.transfer(pointToVerts);
    meshCells_.transfer(cellLabels);
    pointToFace_.transfer(pointToFace);

    tmp<pointField> allPoints
    (
        interpolate
        (
            mesh_.cellCentres(),
            mesh_.points()
        )
    );


    // Assign to MeshedSurface
    faceList allTris(faceLabels.size());
    label verti = 0;
    for (face& allTri : allTris)
    {
        allTri.setSize(3);
        allTri[0] = verts[verti++];
        allTri[1] = verts[verti++];
        allTri[2] = verts[verti++];
    }


1329
    surfZoneList allZones(one(), surfZone("allFaces", allTris.size()));
1330

1331
    MeshStorage::clear();
1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354
    MeshStorage updated
    (
        std::move(allPoints),
        std::move(allTris),
        std::move(allZones)
    );
    MeshStorage::transfer(updated);

    // Now:
    // - generated faces and points are assigned to *this
    // - per point we know:
    //  - pointOnDiag: whether it is on a face-diagonal edge
    //  - pointToFace_: from what pyramid (cell+face) it was produced
    //    (note that the pyramid faces are shared between multiple mesh faces)
    //  - pointToVerts_ : originating mesh vertex or cell centre


    if (debug)
    {
        Pout<< "isoSurfaceTopo : generated " << size() << " faces." << endl;
    }


1355
    if (filter != filterType::NONE)
1356 1357 1358 1359 1360 1361 1362 1363 1364
    {
        // Triangulate outside (filter edges to cell centres and optionally
        // face diagonals)
        DynamicList<label> pointCompactMap(size()); // Back to original point
        DynamicList<label> compactCellIDs(size());  // Per tri the cell
        MeshStorage::operator=
        (
            removeInsidePoints
            (
1365
                (filter == filterType::DIAGCELL ? true : false),
1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387
                *this,
                pointFromDiag,
                pointToFace_,
                startTri,
                pointCompactMap,
                compactCellIDs
            )
        );

        pointToVerts_ = UIndirectList<edge>(pointToVerts_, pointCompactMap)();
        pointToFace_ = UIndirectList<label>(pointToFace_, pointCompactMap)();
        pointFromDiag = UIndirectList<bool>(pointFromDiag, pointCompactMap)();
        meshCells_.transfer(compactCellIDs);

        if (debug)
        {
            Pout<< "isoSurfaceTopo :"
                << " after removing cell centre and face-diag triangles : "
                << size() << endl;
        }


1388
        if (filter == filterType::DIAGCELL)
1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399
        {
            // We remove verts on face diagonals. This is in fact just
            // straightening the edges of the face through the cell. This can
            // close off 'pockets' of triangles and create open or
            // multiply-connected triangles

            // Solved by eroding open-edges
            // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~

            // Mark points on mesh outside. Note that we extend with nCells
            // so we can easily index with pointToVerts_.
1400
            bitSet isBoundaryPoint(mesh.nPoints() + mesh.nCells());
1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411
            for
            (
                label facei = mesh.nInternalFaces();
                facei < mesh.nFaces();
                ++facei
            )
            {
                isBoundaryPoint.set(mesh.faces()[facei]);
            }


1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443
            // Include faces that would be exposed from mesh subset
            if (!ignoreCells_.empty())
            {
                const labelList& faceOwner = mesh_.faceOwner();
                const labelList& faceNeighbour = mesh_.faceNeighbour();

                label nMore = 0;
                for
                (
                    label facei = 0;
                    facei < mesh.nInternalFaces();
                    ++facei
                )
                {
                    int n = 0;
                    if (ignoreCells_.test(faceOwner[facei])) ++n;
                    if (ignoreCells_.test(faceNeighbour[facei])) ++n;

                    // Only one of the cells is being ignored.
                    // That means it is an exposed subMesh face.
                    if (n == 1)
                    {
                        isBoundaryPoint.set(mesh.faces()[facei]);

                        ++nMore;
                    }
                }
            }


            bitSet faceSelection;

1444 1445
            while (true)
            {
1446 1447 1448
                faceSelection.clear();
                faceSelection.resize(this->size());

1449 1450
                const labelList& mp = meshPoints();

1451 1452 1453
                const labelListList& edgeFaces = MeshStorage::edgeFaces();

                forAll(edgeFaces, edgei)
1454
                {
1455 1456
                    const labelList& eFaces = edgeFaces[edgei];
                    if (eFaces.size() == 1)
1457
                    {
1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469
                        // Open edge. Check that vertices do not originate
                        // from a boundary face
                        const edge& e = edges()[edgei];
                        const edge& verts0 = pointToVerts_[mp[e[0]]];
                        const edge& verts1 = pointToVerts_[mp[e[1]]];
                        if
                        (
                            isBoundaryPoint[verts0[0]]
                         && isBoundaryPoint[verts0[1]]
                         && isBoundaryPoint[verts1[0]]
                         && isBoundaryPoint[verts1[1]]
                        )
1470
                        {
1471 1472 1473 1474 1475 1476
                            // Open edge on boundary face. Keep
                        }
                        else
                        {
                            // Open edge. Mark for erosion
                            faceSelection.set(eFaces[0]);
1477 1478 1479 1480 1481 1482 1483
                        }
                    }
                }

                if (debug)
                {
                    Pout<< "isoSurfaceTopo :"
1484
                        << " removing " << faceSelection.count()
1485 1486 1487
                        << " faces since on open edges" << endl;
                }

1488
                if (returnReduce(faceSelection.none(), andOp<bool>()))
1489 1490 1491 1492
                {
                    break;
                }

1493 1494 1495

                // Flip from remove face -> keep face
                faceSelection.flip();
1496 1497 1498 1499 1500 1501 1502

                labelList pointMap;
                labelList faceMap;
                MeshStorage filteredSurf
                (
                    MeshStorage::subsetMesh
                    (
1503
                        faceSelection,
1504 1505 1506 1507 1508 1509 1510 1511 1512 1513 1514 1515 1516 1517 1518 1519 1520
                        pointMap,
                        faceMap
                    )
                );
                MeshStorage::transfer(filteredSurf);

                pointToVerts_ = UIndirectList<edge>(pointToVerts_, pointMap)();
                pointToFace_ = UIndirectList<label>(pointToFace_, pointMap)();
                pointFromDiag = UIndirectList<bool>(pointFromDiag, pointMap)();
                meshCells_ = UIndirectList<label>(meshCells_, faceMap)();
            }
        }
    }
}


// ************************************************************************* //