Skip to content
Snippets Groups Projects
cvMeshBackgroundMesh.C 20.8 KiB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 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 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 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 366 367 368 369 370 371 372 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 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 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
/*---------------------------------------------------------------------------*\
 =========                   |
 \\      /   F ield          | OpenFOAM: The Open Source CFD Toolbox
  \\    /    O peration      |
   \\  /     A nd            | Copyright (C) 2012 OpenFOAM Foundation
    \\/      M anipulation   |
-------------------------------------------------------------------------------
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/>.

Application
    cvMeshBackGroundMesh

Description
    Writes out background mesh as constructed by cvMesh and constructs
    distanceSurface.

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

#include "PatchTools.H"
#include "argList.H"
#include "Time.H"
#include "triSurface.H"
#include "searchableSurfaces.H"
#include "conformationSurfaces.H"
#include "cellSizeControlSurfaces.H"
#include "backgroundMeshDecomposition.H"
#include "cellShape.H"
#include "cellModeller.H"
#include "DynamicField.H"
#include "isoSurfaceCell.H"
#include "vtkSurfaceWriter.H"
#include "syncTools.H"

using namespace Foam;

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

// Tolerance (as fraction of the bounding box). Needs to be fairly lax since
// usually meshes get written with limited precision (6 digits)
static const scalar defaultMergeTol = 1E-6;

// Get merging distance when matching face centres
scalar getMergeDistance
(
    const argList& args,
    const Time& runTime,
    const boundBox& bb
)
{
    scalar mergeTol = defaultMergeTol;
    args.optionReadIfPresent("mergeTol", mergeTol);

    scalar writeTol =
        Foam::pow(scalar(10.0), -scalar(IOstream::defaultPrecision()));

    Info<< "Merge tolerance : " << mergeTol << nl
        << "Write tolerance : " << writeTol << endl;

    if (runTime.writeFormat() == IOstream::ASCII && mergeTol < writeTol)
    {
        FatalErrorIn("getMergeDistance")
            << "Your current settings specify ASCII writing with "
            << IOstream::defaultPrecision() << " digits precision." << endl
            << "Your merging tolerance (" << mergeTol << ") is finer than this."
            << endl
            << "Please change your writeFormat to binary"
            << " or increase the writePrecision" << endl
            << "or adjust the merge tolerance (-mergeTol)."
            << exit(FatalError);
    }

    scalar mergeDist = mergeTol * bb.mag();

    Info<< "Overall meshes bounding box : " << bb << nl
        << "Relative tolerance          : " << mergeTol << nl
        << "Absolute matching distance  : " << mergeDist << nl
        << endl;

    return mergeDist;
}


void printMeshData(const polyMesh& mesh)
{
    // Collect all data on master

    globalIndex globalCells(mesh.nCells());
    labelListList patchNeiProcNo(Pstream::nProcs());
    labelListList patchSize(Pstream::nProcs());
    const labelList& pPatches = mesh.globalData().processorPatches();
    patchNeiProcNo[Pstream::myProcNo()].setSize(pPatches.size());
    patchSize[Pstream::myProcNo()].setSize(pPatches.size());
    forAll(pPatches, i)
    {
        const processorPolyPatch& ppp = refCast<const processorPolyPatch>
        (
            mesh.boundaryMesh()[pPatches[i]]
        );
        patchNeiProcNo[Pstream::myProcNo()][i] = ppp.neighbProcNo();
        patchSize[Pstream::myProcNo()][i] = ppp.size();
    }
    Pstream::gatherList(patchNeiProcNo);
    Pstream::gatherList(patchSize);


    // Print stats

    globalIndex globalBoundaryFaces(mesh.nFaces()-mesh.nInternalFaces());

    label maxProcCells = 0;
    label totProcFaces = 0;
    label maxProcPatches = 0;
    label totProcPatches = 0;
    label maxProcFaces = 0;

    for (label procI = 0; procI < Pstream::nProcs(); procI++)
    {
        Info<< endl
            << "Processor " << procI << nl
            << "    Number of cells = " << globalCells.localSize(procI)
            << endl;

        label nProcFaces = 0;

        const labelList& nei = patchNeiProcNo[procI];

        forAll(patchNeiProcNo[procI], i)
        {
            Info<< "    Number of faces shared with processor "
                << patchNeiProcNo[procI][i] << " = " << patchSize[procI][i]
                << endl;

            nProcFaces += patchSize[procI][i];
        }

        Info<< "    Number of processor patches = " << nei.size() << nl
            << "    Number of processor faces = " << nProcFaces << nl
            << "    Number of boundary faces = "
            << globalBoundaryFaces.localSize(procI) << endl;

        maxProcCells = max(maxProcCells, globalCells.localSize(procI));
        totProcFaces += nProcFaces;
        totProcPatches += nei.size();
        maxProcPatches = max(maxProcPatches, nei.size());
        maxProcFaces = max(maxProcFaces, nProcFaces);
    }

    // Stats

    scalar avgProcCells = scalar(globalCells.size())/Pstream::nProcs();
    scalar avgProcPatches = scalar(totProcPatches)/Pstream::nProcs();
    scalar avgProcFaces = scalar(totProcFaces)/Pstream::nProcs();

    // In case of all faces on one processor. Just to avoid division by 0.
    if (totProcPatches == 0)
    {
        avgProcPatches = 1;
    }
    if (totProcFaces == 0)
    {
        avgProcFaces = 1;
    }

    Info<< nl
        << "Number of processor faces = " << totProcFaces/2 << nl
        << "Max number of cells = " << maxProcCells
        << " (" << 100.0*(maxProcCells-avgProcCells)/avgProcCells
        << "% above average " << avgProcCells << ")" << nl
        << "Max number of processor patches = " << maxProcPatches
        << " (" << 100.0*(maxProcPatches-avgProcPatches)/avgProcPatches
        << "% above average " << avgProcPatches << ")" << nl
        << "Max number of faces between processors = " << maxProcFaces
        << " (" << 100.0*(maxProcFaces-avgProcFaces)/avgProcFaces
        << "% above average " << avgProcFaces << ")" << nl
        << endl;
}


// Return cellID
label cellLabel
(
    const Vector<label>& nCells,
    const label i,
    const label j,
    const label k
)
{
    return i*nCells[1]*nCells[2]+j*nCells[2]+k;
}

label vtxLabel
(
    const Vector<label>& nCells,
    const label i,
    const label j,
    const label k
)
{
    Vector<label> nPoints
    (
        nCells[0]+1,
        nCells[1]+1,
        nCells[2]+1
    );
    return i*nPoints[1]*nPoints[2]+j*nPoints[2]+k;
}


autoPtr<polyMesh> generateHexMesh
(
    const IOobject& io,
    const point& origin,
    const vector& cellSize,
    const Vector<label>& nCells
)
{
    pointField points;
    if (nCells[0]+nCells[1]+nCells[2] > 0)
    {
        points.setSize((nCells[0]+1)*(nCells[1]+1)*(nCells[2]+1));

        // Generate points
        for (label i = 0; i <= nCells[0]; i++)
        {
            for (label j = 0; j <= nCells[1]; j++)
            {
                for (label k = 0; k <= nCells[2]; k++)
                {
                    point pt = origin;
                    pt.x() += i*cellSize[0];
                    pt.y() += j*cellSize[1];
                    pt.z() += k*cellSize[2];
                    points[vtxLabel(nCells, i, j, k)] = pt;
                }
            }
        }
    }


    const cellModel& hex = *(cellModeller::lookup("hex"));
    cellShapeList cellShapes(nCells[0]*nCells[1]*nCells[2]);

    labelList hexPoints(8);
    label cellI = 0;
    for (label i = 0; i < nCells[0]; i++)
    {
        for (label j = 0; j < nCells[1]; j++)
        {
            for (label k = 0; k < nCells[2]; k++)
            {
                hexPoints[0] = vtxLabel(nCells, i,   j,   k);
                hexPoints[1] = vtxLabel(nCells, i+1, j,   k);
                hexPoints[2] = vtxLabel(nCells, i+1, j+1, k);
                hexPoints[3] = vtxLabel(nCells, i,   j+1, k);
                hexPoints[4] = vtxLabel(nCells, i,   j,   k+1);
                hexPoints[5] = vtxLabel(nCells, i+1, j,   k+1);
                hexPoints[6] = vtxLabel(nCells, i+1, j+1, k+1);
                hexPoints[7] = vtxLabel(nCells, i,   j+1, k+1);
                cellShapes[cellI++] = cellShape(hex, hexPoints);
            }
        }
    }

    faceListList boundary(0);
    wordList patchNames(0);
    wordList patchTypes(0);
    word defaultFacesName = "defaultFaces";
    word defaultFacesType = polyPatch::typeName;
    wordList patchPhysicalTypes(0);

    return autoPtr<polyMesh>
    (
        new polyMesh
        (
            io,
            xferMoveTo<pointField>(points),
            cellShapes,
            boundary,
            patchNames,
            patchTypes,
            defaultFacesName,
            defaultFacesType,
            patchPhysicalTypes
        )
    );
}


// Determine for every point a signed distance to the nearest surface
// (outside is positive)
tmp<scalarField> signedDistance
(
    const scalarField& distSqr,
    const pointField& points,
    const searchableSurfaces& geometry,
    const labelList& surfaces
)
{
    tmp<scalarField> tfld(new scalarField(points.size(), Foam::sqr(GREAT)));
    scalarField& fld = tfld();

    // Find nearest
    List<pointIndexHit> nearest;
    labelList nearestSurfaces;
    searchableSurfacesQueries::findNearest
    (
        geometry,
        surfaces,
        points,
        scalarField(points.size(), Foam::sqr(GREAT)),//distSqr
        nearestSurfaces,
        nearest
    );

    // Determine sign of nearest. Sort by surface to do this.
    DynamicField<point> surfPoints(points.size());
    DynamicList<label> surfIndices(points.size());

    forAll(surfaces, surfI)
    {
        // Extract points on this surface
        surfPoints.clear();
        surfIndices.clear();
        forAll(nearestSurfaces, i)
        {
            if (nearestSurfaces[i] == surfI)
            {
                surfPoints.append(points[i]);
                surfIndices.append(i);
            }
        }

        // Calculate sideness of these surface points
        label geomI = surfaces[surfI];
        List<searchableSurface::volumeType> volType;
        geometry[geomI].getVolumeType(surfPoints, volType);

        // Push back to original
        forAll(volType, i)
        {
            label pointI = surfIndices[i];
            scalar dist = mag(points[pointI] - nearest[pointI].hitPoint());

            searchableSurface::volumeType vT = volType[i];

            if (vT == searchableSurface::OUTSIDE)
            {
                fld[pointI] = dist;
            }
            else if (vT == searchableSurface::INSIDE)
            {
                fld[i] = -dist;
            }
            else
            {
                FatalErrorIn("signedDistance()")
                    << "getVolumeType failure, neither INSIDE or OUTSIDE"
                    << exit(FatalError);
            }
        }
    }

    return tfld;
}



// Main program:

int main(int argc, char *argv[])
{
    argList::addNote
    (
        "Generate cvMesh-consistent representation of surfaces"
    );
    argList::addBoolOption
    (
        "writeMesh",
        "write the resulting mesh and distance fields"
    );
    argList::addOption
    (
        "mergeTol",
        "scalar",
        "specify the merge distance relative to the bounding box size "
        "(default 1E-6)"
    );

    #include "setRootCase.H"
    #include "createTime.H"
    runTime.functionObjects().off();

    const bool writeMesh = args.optionFound("writeMesh");

    if (writeMesh)
    {
        Info<< "Writing resulting mesh and cellDistance, pointDistance fields."
            << nl << endl;
    }


    IOdictionary cvMeshDict
    (
        IOobject
        (
            "cvMeshDict",
            runTime.system(),
            runTime,
            IOobject::MUST_READ_IF_MODIFIED,
            IOobject::NO_WRITE
        )
    );

    // Define/load all geometry
    searchableSurfaces allGeometry
    (
        IOobject
        (
            "cvSearchableSurfaces",
            runTime.constant(),
            "triSurface",
            runTime,
            IOobject::MUST_READ,
            IOobject::NO_WRITE
        ),
        cvMeshDict.subDict("geometry")
    );

    Random rndGen(64293*Pstream::myProcNo());

    conformationSurfaces geometryToConformTo
    (
        runTime,
        rndGen,
        allGeometry,
        cvMeshDict.subDict("surfaceConformation")
    );

    cellSizeControlSurfaces cellSizeControl
    (
        allGeometry,
        cvMeshDict.subDict("motionControl")
    );


    // Generate starting block mesh
    vector cellSize;
    {
        const treeBoundBox& bb = geometryToConformTo.globalBounds();

        // Determine the number of cells in each direction.
        const vector span = bb.span();
        vector nScalarCells = span/cellSizeControl.defaultCellSize();

        // Calculate initial cell size to be a little bit smaller than the
        // defaultCellSize to avoid initial refinement triggering.
        Vector<label> nCells = Vector<label>
        (
            label(nScalarCells.x())+2,
            label(nScalarCells.y())+2,
            label(nScalarCells.z())+2
        );
        cellSize = vector
        (
            span[0]/nCells[0],
            span[1]/nCells[1],
            span[2]/nCells[2]
        );

        Info<< "Generating initial hex mesh with" << nl
            << "    bounding box : " << bb << nl
            << "    nCells       : " << nCells << nl
            << "    cellSize     : " << cellSize << nl
            << endl;

        autoPtr<polyMesh> meshPtr
        (
            generateHexMesh
            (
                IOobject
                (
                    polyMesh::defaultRegion,
                    runTime.constant(),
                    runTime
                ),
                bb.min(),
                cellSize,
                (
                    Pstream::master()
                  ? nCells
                  : Vector<label>(0, 0, 0)
                )
            )
        );
        Info<< "Writing initial hex mesh to " << meshPtr().instance() << nl
            << endl;
        meshPtr().write();
    }

    // Distribute the initial mesh
    if (Pstream::parRun())
    {
#       include "createMesh.H"
        Info<< "Loaded mesh:" << endl;
        printMeshData(mesh);

        // Allocate a decomposer
        IOdictionary decompositionDict
        (
            IOobject
            (
                "decomposeParDict",
                runTime.system(),
                mesh,
                IOobject::MUST_READ_IF_MODIFIED,
                IOobject::NO_WRITE
            )
        );

        autoPtr<decompositionMethod> decomposer
        (
            decompositionMethod::New
            (
                decompositionDict
            )
        );

        labelList decomp = decomposer().decompose(mesh, mesh.cellCentres());

        // Global matching tolerance
        const scalar tolDim = getMergeDistance
        (
            args,
            runTime,
            mesh.bounds()
        );

        // Mesh distribution engine
        fvMeshDistribute distributor(mesh, tolDim);

        Info<< "Wanted distribution:"
            << distributor.countCells(decomp) << nl << endl;

        // Do actual sending/receiving of mesh
        autoPtr<mapDistributePolyMesh> map = distributor.distribute(decomp);

        // Print some statistics
        //Info<< "After distribution:" << endl;
        //printMeshData(mesh);

        mesh.setInstance(runTime.constant());
        Info<< "Writing redistributed mesh" << nl << endl;
        mesh.write();
    }


    Info<< "Refining backgroud mesh according to cell size specification" << nl
        << endl;

    backgroundMeshDecomposition backgroundMesh
    (
        1.0,    //spanScale,ratio of poly cell size v.s. hex cell size
        0.0,    //minCellSizeLimit
        0,      //minLevels
        4,      //volRes, check multiple points per cell
        20.0,   //maxCellWeightCoeff
        runTime,
        geometryToConformTo,
        cellSizeControl,
        rndGen,
        cvMeshDict
    );

    if (writeMesh)
    {
        runTime++;
        Info<< "Writing mesh to " << runTime.timeName() << endl;
        backgroundMesh.mesh().write();
    }

    const scalar tolDim = getMergeDistance
    (
        args,
        runTime,
        backgroundMesh.mesh().bounds()
    );


    faceList isoFaces;
    pointField isoPoints;

    {
        // Apply a distanceSurface to it.
        const fvMesh& fvm = backgroundMesh.mesh();

        volScalarField cellDistance
        (
            IOobject
            (
                "cellDistance",
                fvm.time().timeName(),
                fvm.time(),
                IOobject::NO_READ,
                IOobject::NO_WRITE,
                false
            ),
            fvm,
            dimensionedScalar("zero", dimLength, 0)
        );

        const searchableSurfaces& geometry = geometryToConformTo.geometry();
        const labelList& surfaces = geometryToConformTo.surfaces();


        // Get maximum search size per cell
        scalarField distSqr(cellDistance.size());

        const labelList& cellLevel = backgroundMesh.cellLevel();
        forAll(cellLevel, cellI)
        {
            // The largest edge of the cell will always be less than the
            // span of the bounding box of the cell.
            distSqr[cellI] = magSqr(cellSize)/pow(2, cellLevel[cellI]);
        }

        {
            // Internal field
            cellDistance.internalField() = signedDistance
            (
                distSqr,
                fvm.C(),
                geometry,
                surfaces
            );
            // Patch fields
            forAll(fvm.C().boundaryField(), patchI)
            {
                const pointField& cc = fvm.C().boundaryField()[patchI];
                fvPatchScalarField& fld = cellDistance.boundaryField()[patchI];
                scalarField patchDistSqr
                (
                    fld.patch().patchInternalField(distSqr)
                );
                fld = signedDistance(patchDistSqr, cc, geometry, surfaces);
            }

            // On processor patches the fvm.C() will already be the cell centre
            // on the opposite side so no need to swap cellDistance.

            if (writeMesh)
            {
                cellDistance.write();
            }
        }


        // Distance to points
        pointScalarField pointDistance
        (
            IOobject
            (
                "pointDistance",
                fvm.time().timeName(),
                fvm.time(),
                IOobject::NO_READ,
                IOobject::NO_WRITE,
                false
            ),
            pointMesh::New(fvm),
            dimensionedScalar("zero", dimLength, 0)
        );
        {
            scalarField pointDistSqr(fvm.nPoints(), -sqr(GREAT));
            for (label faceI = 0; faceI < fvm.nInternalFaces(); faceI++)
            {
                label own = fvm.faceOwner()[faceI];
                label ownDistSqr = distSqr[own];

                const face& f = fvm.faces()[faceI];
                forAll(f, fp)
                {
                    pointDistSqr[f[fp]] = max(pointDistSqr[f[fp]], ownDistSqr);
                }
            }
            syncTools::syncPointList
            (
                fvm,
                pointDistSqr,
                maxEqOp<scalar>(),
                -sqr(GREAT)             // null value
            );

            pointDistance.internalField() = signedDistance
            (
                pointDistSqr,
                fvm.points(),
                geometry,
                surfaces
            );

            if (writeMesh)
            {
                pointDistance.write();
            }
        }

        isoSurfaceCell iso
        (
            fvm,
            cellDistance,
            pointDistance,
            0,      //distance,
            false   //regularise
        );

        isoFaces.setSize(iso.size());
        forAll(isoFaces, i)
        {
            isoFaces[i] = iso[i].triFaceFace();
        }
        isoPoints = iso.points();
    }


    pointField mergedPoints;
    faceList mergedFaces;
    labelList pointMergeMap;
    PatchTools::gatherAndMerge
    (
        tolDim,
        primitivePatch
        (
            SubList<face>(isoFaces, isoFaces.size()),
            isoPoints
        ),
        mergedPoints,
        mergedFaces,
        pointMergeMap
    );

    vtkSurfaceWriter writer;
    writer.write
    (
        runTime.path(),
        "iso",
        mergedPoints,
        mergedFaces
    );

    Info<< "End\n" << endl;

    return 0;
}


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