multiphaseSystem.C 17.7 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2013-2017 OpenFOAM Foundation
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
     \\/     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/>.

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

#include "multiphaseSystem.H"
#include "alphaContactAngleFvPatchScalarField.H"

#include "MULES.H"
#include "subCycle.H"

#include "fvcDdt.H"
#include "fvcDiv.H"
#include "fvcSnGrad.H"
#include "fvcFlux.H"
#include "fvcMeshPhi.H"
#include "fvcSup.H"

#include "fvmDdt.H"
#include "fvmLaplacian.H"
#include "fvmSup.H"

43 44
#include "unitConversion.H"

45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60
// * * * * * * * * * * * * * * * Static Member Data  * * * * * * * * * * * * //

namespace Foam
{
    defineTypeNameAndDebug(multiphaseSystem, 0);
    defineRunTimeSelectionTable(multiphaseSystem, dictionary);
}


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

void Foam::multiphaseSystem::calcAlphas()
{
    scalar level = 0.0;
    alphas_ == 0.0;

61
    for (const phaseModel& phase : phases())
62
    {
63
        alphas_ += level * phase;
64 65 66 67 68 69 70
        level += 1.0;
    }
}


void Foam::multiphaseSystem::solveAlphas()
{
71 72
    bool LTS = fv::localEulerDdt::enabled(mesh_);

73
    for (phaseModel& phase : phases())
74
    {
75
        phase.correctBoundaryConditions();
76 77
    }

78
    PtrList<surfaceScalarField> alphaPhiCorrs(phases().size());
79 80 81

    label phasei = 0;
    for (phaseModel& phase : phases())
82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101
    {
        volScalarField& alpha1 = phase;

        alphaPhiCorrs.set
        (
            phasei,
            new surfaceScalarField
            (
                "phi" + alpha1.name() + "Corr",
                fvc::flux
                (
                    phi_,
                    phase,
                    "div(phi," + alpha1.name() + ')'
                )
            )
        );

        surfaceScalarField& alphaPhiCorr = alphaPhiCorrs[phasei];

102
        for (const phaseModel& phase2 : phases())
103
        {
104
            const volScalarField& alpha2 = phase2;
105 106 107 108 109

            if (&phase2 == &phase) continue;

            surfaceScalarField phir(phase.phi() - phase2.phi());

110 111
            const auto cAlpha =
                cAlphas_.cfind(phasePairKey(phase.name(), phase2.name()));
112

113
            if (cAlpha.found())
114 115 116 117 118 119 120 121 122
            {
                surfaceScalarField phic
                (
                    (mag(phi_) + mag(phir))/mesh_.magSf()
                );

                phir += min(cAlpha()*phic, max(phic))*nHatf(phase, phase2);
            }

123
            const word phirScheme
124 125 126 127 128 129 130 131 132 133 134 135
            (
                "div(phir," + alpha2.name() + ',' + alpha1.name() + ')'
            );

            alphaPhiCorr += fvc::flux
            (
                -fvc::flux(-phir, phase2, phirScheme),
                phase,
                phirScheme
            );
        }

136
        phase.correctInflowOutflow(alphaPhiCorr);
137

138
        if (LTS)
139 140 141
        {
            MULES::limit
            (
142
                fv::localEulerDdt::localRDeltaT(mesh_),
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
                geometricOneField(),
                phase,
                phi_,
                alphaPhiCorr,
                zeroField(),
                zeroField(),
                phase.alphaMax(),
                0,
                true
            );
        }
        else
        {
            const scalar rDeltaT = 1.0/mesh_.time().deltaTValue();

            MULES::limit
            (
                rDeltaT,
                geometricOneField(),
                phase,
                phi_,
                alphaPhiCorr,
                zeroField(),
                zeroField(),
                phase.alphaMax(),
                0,
                true
            );
        }
172 173

        ++phasei;
174 175 176 177 178 179 180 181 182 183 184 185 186
    }

    MULES::limitSum(alphaPhiCorrs);

    volScalarField sumAlpha
    (
        IOobject
        (
            "sumAlpha",
            mesh_.time().timeName(),
            mesh_
        ),
        mesh_,
187
        dimensionedScalar(dimless, Zero)
188 189 190 191 192
    );


    volScalarField divU(fvc::div(fvc::absolute(phi_, phases().first().U())));

193
    forAll(phases(), phasei)
194
    {
195
        phaseModel& phase = phases()[phasei];
196 197
        volScalarField& alpha = phase;

198 199 200
        surfaceScalarField& alphaPhi = alphaPhiCorrs[phasei];
        alphaPhi += upwind<scalar>(mesh_, phi_).flux(phase);
        phase.correctInflowOutflow(alphaPhi);
201

202
        volScalarField::Internal Sp
203 204 205 206 207 208 209 210
        (
            IOobject
            (
                "Sp",
                mesh_.time().timeName(),
                mesh_
            ),
            mesh_,
211
            dimensionedScalar(divU.dimensions(), Zero)
212 213
        );

214
        volScalarField::Internal Su
215 216 217 218 219 220 221 222 223 224 225 226
        (
            IOobject
            (
                "Su",
                mesh_.time().timeName(),
                mesh_
            ),
            // Divergence term is handled explicitly to be
            // consistent with the explicit transport solution
            divU*min(alpha, scalar(1))
        );

227
        if (phase.divU().valid())
228
        {
229
            const scalarField& dgdt = phase.divU()();
230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246

            forAll(dgdt, celli)
            {
                if (dgdt[celli] > 0.0)
                {
                    Sp[celli] -= dgdt[celli];
                    Su[celli] += dgdt[celli];
                }
                else if (dgdt[celli] < 0.0)
                {
                    Sp[celli] +=
                        dgdt[celli]
                       *(1.0 - alpha[celli])/max(alpha[celli], 1e-4);
                }
            }
        }

247
        for (const phaseModel& phase2 : phases())
248 249 250 251 252
        {
            const volScalarField& alpha2 = phase2;

            if (&phase2 == &phase) continue;

253
            if (phase2.divU().valid())
254
            {
255
                const scalarField& dgdt2 = phase2.divU()();
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

                forAll(dgdt2, celli)
                {
                    if (dgdt2[celli] < 0.0)
                    {
                        Sp[celli] +=
                            dgdt2[celli]
                           *(1.0 - alpha2[celli])/max(alpha2[celli], 1e-4);

                        Su[celli] -=
                            dgdt2[celli]
                           *alpha[celli]/max(alpha2[celli], 1e-4);
                    }
                    else if (dgdt2[celli] > 0.0)
                    {
                        Sp[celli] -= dgdt2[celli];
                    }
                }
            }
        }

        MULES::explicitSolve
        (
            geometricOneField(),
            alpha,
281
            alphaPhi,
282 283 284 285
            Sp,
            Su
        );

286
        phase.alphaPhi() = alphaPhi;
287 288 289 290 291 292 293 294 295 296 297 298 299 300 301

        Info<< phase.name() << " volume fraction, min, max = "
            << phase.weightedAverage(mesh_.V()).value()
            << ' ' << min(phase).value()
            << ' ' << max(phase).value()
            << endl;

        sumAlpha += phase;
    }

    Info<< "Phase-sum volume fraction, min, max = "
        << sumAlpha.weightedAverage(mesh_.V()).value()
        << ' ' << min(sumAlpha).value()
        << ' ' << max(sumAlpha).value()
        << endl;
302 303 304

    // Correct the sum of the phase-fractions to avoid 'drift'
    volScalarField sumCorr(1.0 - sumAlpha);
305
    for (phaseModel& phase : phases())
306
    {
307
        volScalarField& alpha = phase;
308 309
        alpha += alpha*sumCorr;
    }
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
}


Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseSystem::nHatfv
(
    const volScalarField& alpha1,
    const volScalarField& alpha2
) const
{
    /*
    // Cell gradient of alpha
    volVectorField gradAlpha =
        alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);

    // Interpolated face-gradient of alpha
    surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
    */

    surfaceVectorField gradAlphaf
    (
        fvc::interpolate(alpha2)*fvc::interpolate(fvc::grad(alpha1))
      - fvc::interpolate(alpha1)*fvc::interpolate(fvc::grad(alpha2))
    );

    // Face unit interface normal
    return gradAlphaf/(mag(gradAlphaf) + deltaN_);
}


Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseSystem::nHatf
(
    const volScalarField& alpha1,
    const volScalarField& alpha2
) const
{
    // Face unit interface normal flux
    return nHatfv(alpha1, alpha2) & mesh_.Sf();
}


// Correction for the boundary condition on the unit normal nHat on
// walls to produce the correct contact angle.

// The dynamic contact angle is calculated from the component of the
// velocity on the direction of the interface, parallel to the wall.

void Foam::multiphaseSystem::correctContactAngle
(
    const phaseModel& phase1,
    const phaseModel& phase2,
360
    surfaceVectorField::Boundary& nHatb
361 362
) const
{
363
    const volScalarField::Boundary& gbf
364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382
        = phase1.boundaryField();

    const fvBoundaryMesh& boundary = mesh_.boundary();

    forAll(boundary, patchi)
    {
        if (isA<alphaContactAngleFvPatchScalarField>(gbf[patchi]))
        {
            const alphaContactAngleFvPatchScalarField& acap =
                refCast<const alphaContactAngleFvPatchScalarField>(gbf[patchi]);

            vectorField& nHatPatch = nHatb[patchi];

            vectorField AfHatPatch
            (
                mesh_.Sf().boundaryField()[patchi]
               /mesh_.magSf().boundaryField()[patchi]
            );

383 384 385
            const auto tp =
                acap.thetaProps()
                    .cfind(phasePairKey(phase1.name(), phase2.name()));
386

387
            if (!tp.found())
388
            {
389 390
                FatalErrorInFunction
                    << "Cannot find interface "
391 392 393 394 395 396
                    << phasePairKey(phase1.name(), phase2.name())
                    << "\n    in table of theta properties for patch "
                    << acap.patch().name()
                    << exit(FatalError);
            }

397
            const bool matched = (tp.key().first() == phase1.name());
398

399
            const scalar theta0 = degToRad(tp().theta0(matched));
400 401 402 403 404 405 406
            scalarField theta(boundary[patchi].size(), theta0);

            scalar uTheta = tp().uTheta();

            // Calculate the dynamic contact angle if required
            if (uTheta > SMALL)
            {
407 408
                const scalar thetaA = degToRad(tp().thetaA(matched));
                const scalar thetaR = degToRad(tp().thetaR(matched));
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

                // Calculated the component of the velocity parallel to the wall
                vectorField Uwall
                (
                    phase1.U()().boundaryField()[patchi].patchInternalField()
                  - phase1.U()().boundaryField()[patchi]
                );
                Uwall -= (AfHatPatch & Uwall)*AfHatPatch;

                // Find the direction of the interface parallel to the wall
                vectorField nWall
                (
                    nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
                );

                // Normalise nWall
                nWall /= (mag(nWall) + SMALL);

                // Calculate Uwall resolved normal to the interface parallel to
                // the interface
                scalarField uwall(nWall & Uwall);

                theta += (thetaA - thetaR)*tanh(uwall/uTheta);
            }


            // Reset nHatPatch to correspond to the contact angle

            scalarField a12(nHatPatch & AfHatPatch);

            scalarField b1(cos(theta));

            scalarField b2(nHatPatch.size());

            forAll(b2, facei)
            {
                b2[facei] = cos(acos(a12[facei]) - theta[facei]);
            }

            scalarField det(1.0 - a12*a12);

            scalarField a((b1 - a12*b2)/det);
            scalarField b((b2 - a12*b1)/det);

            nHatPatch = a*AfHatPatch + b*nHatPatch;

            nHatPatch /= (mag(nHatPatch) + deltaN_.value());
        }
    }
}


Foam::tmp<Foam::volScalarField> Foam::multiphaseSystem::K
(
    const phaseModel& phase1,
    const phaseModel& phase2
) const
{
    tmp<surfaceVectorField> tnHatfv = nHatfv(phase1, phase2);

469
    correctContactAngle(phase1, phase2, tnHatfv.ref().boundaryFieldRef());
470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489

    // Simple expression for curvature
    return -fvc::div(tnHatfv & mesh_.Sf());
}


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

Foam::multiphaseSystem::multiphaseSystem
(
    const fvMesh& mesh
)
:
    phaseSystem(mesh),

    alphas_
    (
        IOobject
        (
            "alphas",
490
            mesh_.time().timeName(),
491 492 493 494 495
            mesh,
            IOobject::NO_READ,
            IOobject::AUTO_WRITE
        ),
        mesh,
496
        dimensionedScalar(dimless, Zero)
497 498 499 500 501 502 503
    ),

    cAlphas_(lookup("interfaceCompression")),

    deltaN_
    (
        "deltaN",
504
        1e-8/cbrt(average(mesh_.V()))
505 506
    )
{
507
    for (const phaseModel& phase : phases())
508
    {
509 510
        const volScalarField& alpha = phase;
        mesh_.setFluxRequired(alpha.name());
511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532
    }
}


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

Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseSystem::surfaceTension
(
    const phaseModel& phase1
) const
{
    tmp<surfaceScalarField> tSurfaceTension
    (
        new surfaceScalarField
        (
            IOobject
            (
                "surfaceTension",
                mesh_.time().timeName(),
                mesh_
            ),
            mesh_,
533
            dimensionedScalar(dimensionSet(1, -2, -2, 0, 0), Zero)
534 535 536
        )
    );

537
    tSurfaceTension.ref().setOriented();
538

539
    for (const phaseModel& phase2 : phases())
540 541 542 543 544
    {
        if (&phase2 != &phase1)
        {
            phasePairKey key12(phase1.name(), phase2.name());

545
            const auto cAlpha = cAlphas_.cfind(key12);
546

547
            if (cAlpha.found())
548
            {
549
                tSurfaceTension.ref() +=
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
                    fvc::interpolate(sigma(key12)*K(phase1, phase2))
                   *(
                        fvc::interpolate(phase2)*fvc::snGrad(phase1)
                      - fvc::interpolate(phase1)*fvc::snGrad(phase2)
                    );
            }
        }
    }

    return tSurfaceTension;
}


Foam::tmp<Foam::volScalarField>
Foam::multiphaseSystem::nearInterface() const
{
    tmp<volScalarField> tnearInt
    (
        new volScalarField
        (
            IOobject
            (
                "nearInterface",
                mesh_.time().timeName(),
                mesh_
            ),
            mesh_,
577
            dimensionedScalar(dimless, Zero)
578 579 580
        )
    );

581
    for (const phaseModel& phase : phases())
582
    {
583
        tnearInt.ref() = max
584 585
        (
            tnearInt(),
586
            pos0(phase - 0.01)*pos0(0.99 - phase)
587
        );
588 589 590 591 592 593 594 595
    }

    return tnearInt;
}


void Foam::multiphaseSystem::solve()
{
596
    const Time& runTime = mesh_.time();
597 598

    const dictionary& alphaControls = mesh_.solverDict("alpha");
599
    label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
600

601
    bool LTS = fv::localEulerDdt::enabled(mesh_);
602 603 604 605 606 607 608 609

    if (nAlphaSubCycles > 1)
    {
        tmp<volScalarField> trSubDeltaT;

        if (LTS)
        {
            trSubDeltaT =
610
                fv::localEulerDdt::localRSubDeltaT(mesh_, nAlphaSubCycles);
611 612 613
        }

        PtrList<volScalarField> alpha0s(phases().size());
614
        PtrList<surfaceScalarField> alphaPhiSums(phases().size());
615

616 617
        label phasei = 0;
        for (phaseModel& phase : phases())
618 619 620 621 622 623 624 625 626
        {
            volScalarField& alpha = phase;

            alpha0s.set
            (
                phasei,
                new volScalarField(alpha.oldTime())
            );

627
            alphaPhiSums.set
628 629 630 631 632 633 634 635 636 637 638
            (
                phasei,
                new surfaceScalarField
                (
                    IOobject
                    (
                        "phiSum" + alpha.name(),
                        runTime.timeName(),
                        mesh_
                    ),
                    mesh_,
639
                    dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero)
640 641
                )
            );
642 643

            ++phasei;
644 645 646 647 648 649 650 651 652 653 654 655 656 657
        }

        for
        (
            subCycleTime alphaSubCycle
            (
                const_cast<Time&>(runTime),
                nAlphaSubCycles
            );
            !(++alphaSubCycle).end();
        )
        {
            solveAlphas();

658 659
            phasei = 0;
            for (const phaseModel& phase : phases())
660
            {
661 662 663
                alphaPhiSums[phasei] += phase.alphaPhi();

                ++phasei;
664 665 666
            }
        }

667 668
        phasei = 0;
        for (phaseModel& phase : phases())
669 670 671
        {
            volScalarField& alpha = phase;

672
            phase.alphaPhi() = alphaPhiSums[phasei]/nAlphaSubCycles;
673 674 675 676 677 678 679 680

            // Correct the time index of the field
            // to correspond to the global time
            alpha.timeIndex() = runTime.timeIndex();

            // Reset the old-time field value
            alpha.oldTime() = alpha0s[phasei];
            alpha.oldTime().timeIndex() = runTime.timeIndex();
681 682

            ++phasei;
683 684 685 686 687 688 689
        }
    }
    else
    {
        solveAlphas();
    }

690
    for (phaseModel& phase : phases())
691 692
    {
        phase.alphaRhoPhi() = fvc::interpolate(phase.rho())*phase.alphaPhi();
693 694

        // Ensure the phase-fractions are bounded
695
        phase.clip(0, 1);
696 697 698 699 700 701 702
    }

    calcAlphas();
}


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