multiphaseSystem.C 18.7 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
OpenFOAM bot's avatar
OpenFOAM bot committed
5
    \\  /    A nd           | www.openfoam.com
6
     \\/     M anipulation  |
OpenFOAM bot's avatar
OpenFOAM bot committed
7
-------------------------------------------------------------------------------
OpenFOAM bot's avatar
OpenFOAM bot committed
8
    Copyright (C) 2013-2019 OpenFOAM Foundation
9
    Copyright (C) 2020 OpenCFD Ltd.
10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33
-------------------------------------------------------------------------------
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"
34
#include "UniformField.H"
35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54

#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"

// * * * * * * * * * * * * * * * Static Member Data  * * * * * * * * * * * * //

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

55 56 57
const Foam::scalar Foam::multiphaseSystem::convertToRad =
    Foam::constant::mathematical::pi/180.0;

58 59 60 61 62 63 64 65

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

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

sergio's avatar
sergio committed
66
    forAll(phases(), i)
67
    {
sergio's avatar
sergio committed
68
        alphas_ += level*phases()[i];
69 70 71 72 73 74 75
        level += 1.0;
    }
}


void Foam::multiphaseSystem::solveAlphas()
{
sergio's avatar
sergio committed
76
    forAll(phases(), phasei)
77
    {
sergio's avatar
sergio committed
78
        phases()[phasei].correctBoundaryConditions();
79 80
    }

81 82 83 84 85 86 87 88 89 90 91 92
    // Calculate the void fraction
    volScalarField alphaVoid
    (
        IOobject
        (
            "alphaVoid",
            mesh_.time().timeName(),
            mesh_
        ),
        mesh_,
        dimensionedScalar("one", dimless, 1)
    );
sergio's avatar
sergio committed
93
    forAll(stationaryPhases(), stationaryPhasei)
94
    {
sergio's avatar
sergio committed
95
        alphaVoid -= stationaryPhases()[stationaryPhasei];
96
    }
97

98 99 100
    // Generate face-alphas
    PtrList<surfaceScalarField> alphafs(phases().size());
    forAll(phases(), phasei)
101
    {
102 103 104 105 106 107 108 109 110 111 112
        phaseModel& phase = phases()[phasei];
        alphafs.set
        (
            phasei,
            new surfaceScalarField
            (
                IOobject::groupName("alphaf", phase.name()),
                upwind<scalar>(mesh_, phi_).interpolate(phase)
            )
        );
    }
113

114 115
    // Create correction fluxes
    PtrList<surfaceScalarField> alphaPhiCorrs(phases().size());
sergio's avatar
sergio committed
116
    forAll(stationaryPhases(), stationaryPhasei)
117
    {
sergio's avatar
sergio committed
118 119
        phaseModel& phase = stationaryPhases()[stationaryPhasei];

120 121
        alphaPhiCorrs.set
        (
122
            phase.index(),
123 124
            new surfaceScalarField
            (
125 126 127 128 129
                IOobject::groupName("alphaPhiCorr", phase.name()),
              - upwind<scalar>(mesh_, phi_).flux(phase)
            )
        );
    }
sergio's avatar
sergio committed
130
    forAll(movingPhases(), movingPhasei)
131
    {
sergio's avatar
sergio committed
132 133
        phaseModel& phase = movingPhases()[movingPhasei];
        volScalarField& alpha = phase;
134 135 136 137 138 139 140 141

        alphaPhiCorrs.set
        (
            phase.index(),
            new surfaceScalarField
            (
                IOobject::groupName("alphaPhiCorr", phase.name()),
                fvc::flux(phi_, phase, "div(phi," + alpha.name() + ')')
142 143 144
            )
        );

145
        surfaceScalarField& alphaPhiCorr = alphaPhiCorrs[phase.index()];
146

sergio's avatar
sergio committed
147
        forAll(phases(), phasei)
148
        {
sergio's avatar
sergio committed
149 150
            phaseModel& phase2 = phases()[phasei];
            volScalarField& alpha2 = phase2;
151 152 153 154 155

            if (&phase2 == &phase) continue;

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

sergio's avatar
sergio committed
156 157 158 159
            cAlphaTable::const_iterator cAlpha
            (
                cAlphas_.find(phasePairKey(phase.name(), phase2.name()))
            );
160

sergio's avatar
sergio committed
161
            if (cAlpha != cAlphas_.end())
162 163 164 165 166 167 168 169 170
            {
                surfaceScalarField phic
                (
                    (mag(phi_) + mag(phir))/mesh_.magSf()
                );

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

sergio's avatar
sergio committed
171
            word phirScheme
172
            (
173
                "div(phir," + alpha2.name() + ',' + alpha.name() + ')'
174 175 176 177 178 179 180 181 182 183
            );

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

184
        phase.correctInflowOutflow(alphaPhiCorr);
185

186
        MULES::limit
187
        (
188 189 190 191 192 193 194 195 196 197 198
            geometricOneField(),
            phase,
            phi_,
            alphaPhiCorr,
            zeroField(),
            zeroField(),
            min(alphaVoid.primitiveField(), phase.alphaMax())(),
            zeroField(),
            true
        );
    }
199

200 201
    // Limit the flux sums, fixing those of the stationary phases
    labelHashSet fixedAlphaPhiCorrs;
sergio's avatar
sergio committed
202
    forAll(stationaryPhases(), stationaryPhasei)
203
    {
sergio's avatar
sergio committed
204
        fixedAlphaPhiCorrs.insert(stationaryPhases()[stationaryPhasei].index());
205 206
    }
    MULES::limitSum(alphafs, alphaPhiCorrs, fixedAlphaPhiCorrs);
207

208
    // Solve for the moving phase alphas
sergio's avatar
sergio committed
209
    forAll(movingPhases(), movingPhasei)
210
    {
sergio's avatar
sergio committed
211
        phaseModel& phase = movingPhases()[movingPhasei];
212 213
        volScalarField& alpha = phase;

214
        surfaceScalarField& alphaPhi = alphaPhiCorrs[phase.index()];
215 216
        alphaPhi += upwind<scalar>(mesh_, phi_).flux(phase);
        phase.correctInflowOutflow(alphaPhi);
217

218
        volScalarField::Internal Sp
219 220 221 222 223 224 225 226
        (
            IOobject
            (
                "Sp",
                mesh_.time().timeName(),
                mesh_
            ),
            mesh_,
sergio's avatar
sergio committed
227
            dimensionedScalar("zero", dimless/dimTime, 0)
228 229
        );

230
        volScalarField::Internal Su
231
        (
232 233
            "Su",
            min(alpha, scalar(1))*fvc::div(fvc::absolute(phi_, phase.U()))
234 235
        );

236
        if (phase.divU().valid())
237
        {
238
            const scalarField& dgdt = phase.divU()();
239 240 241 242 243 244 245 246 247 248 249 250

            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]
251
                       *(1 - alpha[celli])/max(alpha[celli], 1e-4);
252 253 254 255
                }
            }
        }

sergio's avatar
sergio committed
256
        forAll(phases(), phasej)
257
        {
sergio's avatar
sergio committed
258
            const phaseModel& phase2 = phases()[phasej];
259 260 261 262
            const volScalarField& alpha2 = phase2;

            if (&phase2 == &phase) continue;

263
            if (phase2.divU().valid())
264
            {
265
                const scalarField& dgdt2 = phase2.divU()();
266 267 268 269 270 271 272

                forAll(dgdt2, celli)
                {
                    if (dgdt2[celli] < 0.0)
                    {
                        Sp[celli] +=
                            dgdt2[celli]
273
                           *(1 - alpha2[celli])/max(alpha2[celli], 1e-4);
274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290

                        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,
291
            alphaPhi,
292 293 294 295
            Sp,
            Su
        );

296 297
        phase.alphaPhiRef() = alphaPhi;
    }
298

299
    // Report the phase fractions and the phase fraction sum
sergio's avatar
sergio committed
300
    forAll(phases(), phasei)
301
    {
sergio's avatar
sergio committed
302 303
        phaseModel& phase = phases()[phasei];

304
        Info<< phase.name() << " fraction, min, max = "
305 306 307 308
            << phase.weightedAverage(mesh_.V()).value()
            << ' ' << min(phase).value()
            << ' ' << max(phase).value()
            << endl;
309
    }
310

311 312 313 314 315 316 317 318 319
    volScalarField sumAlphaMoving
    (
        IOobject
        (
            "sumAlphaMoving",
            mesh_.time().timeName(),
            mesh_
        ),
        mesh_,
sergio's avatar
sergio committed
320
        dimensionedScalar("zero", dimless, 0)
321
    );
sergio's avatar
sergio committed
322
    forAll(movingPhases(), movingPhasei)
323
    {
sergio's avatar
sergio committed
324
        sumAlphaMoving += movingPhases()[movingPhasei];
325 326 327
    }

    Info<< "Phase-sum volume fraction, min, max = "
328 329 330
        << (sumAlphaMoving + 1 - alphaVoid)().weightedAverage(mesh_.V()).value()
        << ' ' << min(sumAlphaMoving + 1 - alphaVoid).value()
        << ' ' << max(sumAlphaMoving + 1 - alphaVoid).value()
331
        << endl;
332

333
    // Correct the sum of the phase fractions to avoid drift
sergio's avatar
sergio committed
334
    forAll(movingPhases(), movingPhasei)
335
    {
sergio's avatar
sergio committed
336
        movingPhases()[movingPhasei] *= alphaVoid/sumAlphaMoving;
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
}


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,
388
    surfaceVectorField::Boundary& nHatb
389 390
) const
{
391
    const volScalarField::Boundary& gbf
392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410
        = 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]
            );

sergio's avatar
sergio committed
411 412 413 414
            alphaContactAngleFvPatchScalarField::thetaPropsTable::
                const_iterator tp =
                    acap.thetaProps()
                   .find(phasePairKey(phase1.name(), phase2.name()));
415

sergio's avatar
sergio committed
416
            if (tp == acap.thetaProps().end())
417
            {
418 419
                FatalErrorInFunction
                    << "Cannot find interface "
420 421 422 423 424 425
                    << phasePairKey(phase1.name(), phase2.name())
                    << "\n    in table of theta properties for patch "
                    << acap.patch().name()
                    << exit(FatalError);
            }

sergio's avatar
sergio committed
426
            bool matched = (tp.key().first() == phase1.name());
427

428
            scalar theta0 = convertToRad*tp().theta0(matched);
429 430 431 432 433 434 435
            scalarField theta(boundary[patchi].size(), theta0);

            scalar uTheta = tp().uTheta();

            // Calculate the dynamic contact angle if required
            if (uTheta > SMALL)
            {
436 437
                scalar thetaA = convertToRad*tp().thetaA(matched);
                scalar thetaR = convertToRad*tp().thetaR(matched);
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

                // 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]);
            }

477
            scalarField det(1 - a12*a12);
478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497

            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);

498
    correctContactAngle(phase1, phase2, tnHatfv.ref().boundaryFieldRef());
499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518

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


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

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

    alphas_
    (
        IOobject
        (
            "alphas",
519
            mesh_.time().timeName(),
520 521 522 523 524
            mesh,
            IOobject::NO_READ,
            IOobject::AUTO_WRITE
        ),
        mesh,
sergio's avatar
sergio committed
525
        dimensionedScalar("zero", dimless, 0)
526 527 528 529 530 531 532
    ),

    cAlphas_(lookup("interfaceCompression")),

    deltaN_
    (
        "deltaN",
sergio's avatar
sergio committed
533
        1e-8/pow(average(mesh_.V()), 1.0/3.0)
534 535
    )
{
sergio's avatar
sergio committed
536
    forAll(phases(), phasei)
537
    {
sergio's avatar
sergio committed
538 539
        volScalarField& alphai = phases()[phasei];
        mesh_.setFluxRequired(alphai.name());
540 541 542 543
    }
}


sergio's avatar
sergio committed
544 545 546 547 548 549
// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //

Foam::multiphaseSystem::~multiphaseSystem()
{}


550 551 552 553 554 555 556 557 558
// * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * * //

Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseSystem::surfaceTension
(
    const phaseModel& phase1
) const
{
    tmp<surfaceScalarField> tSurfaceTension
    (
559
        surfaceScalarField::New
560
        (
561
            "surfaceTension",
562
            mesh_,
sergio's avatar
sergio committed
563
            dimensionedScalar("zero", dimensionSet(1, -2, -2, 0, 0), 0)
564 565 566
        )
    );

567
    tSurfaceTension.ref().setOriented();
568

sergio's avatar
sergio committed
569
    forAll(phases(), phasej)
570
    {
sergio's avatar
sergio committed
571 572
        const phaseModel& phase2 = phases()[phasej];

573 574 575 576
        if (&phase2 != &phase1)
        {
            phasePairKey key12(phase1.name(), phase2.name());

sergio's avatar
sergio committed
577
            cAlphaTable::const_iterator cAlpha(cAlphas_.find(key12));
578

sergio's avatar
sergio committed
579
            if (cAlpha != cAlphas_.end())
580
            {
581
                tSurfaceTension.ref() +=
582 583 584 585 586 587 588 589
                    fvc::interpolate(sigma(key12)*K(phase1, phase2))
                   *(
                        fvc::interpolate(phase2)*fvc::snGrad(phase1)
                      - fvc::interpolate(phase1)*fvc::snGrad(phase2)
                    );
            }
        }
    }
sergio's avatar
sergio committed
590

591
    tSurfaceTension->setOriented();
592 593 594 595 596 597 598 599 600 601

    return tSurfaceTension;
}


Foam::tmp<Foam::volScalarField>
Foam::multiphaseSystem::nearInterface() const
{
    tmp<volScalarField> tnearInt
    (
602
        volScalarField::New
603
        (
604
            "nearInterface",
605
            mesh_,
sergio's avatar
sergio committed
606
            dimensionedScalar("zero", dimless, 0)
607 608 609
        )
    );

sergio's avatar
sergio committed
610
    forAll(phases(), phasei)
611
    {
612
        tnearInt.ref() = max
613 614
        (
            tnearInt(),
sergio's avatar
sergio committed
615
            pos0(phases()[phasei] - 0.01)*pos0(0.99 - phases()[phasei])
616
        );
617 618 619 620 621 622 623 624
    }

    return tnearInt;
}


void Foam::multiphaseSystem::solve()
{
625
    const Time& runTime = mesh_.time();
626 627

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

630
    bool LTS = fv::localEulerDdt::enabled(mesh_);
631 632 633 634 635 636 637 638

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

        if (LTS)
        {
            trSubDeltaT =
639
                fv::localEulerDdt::localRSubDeltaT(mesh_, nAlphaSubCycles);
640 641
        }

642
        List<volScalarField*> alphaPtrs(phases().size());
643
        PtrList<surfaceScalarField> alphaPhiSums(phases().size());
644

sergio's avatar
sergio committed
645
        forAll(phases(), phasei)
646
        {
sergio's avatar
sergio committed
647
            phaseModel& phase = phases()[phasei];
648 649
            volScalarField& alpha = phase;

650
            alphaPtrs[phasei] = &alpha;
651

652
            alphaPhiSums.set
653 654 655 656 657 658 659 660 661 662 663
            (
                phasei,
                new surfaceScalarField
                (
                    IOobject
                    (
                        "phiSum" + alpha.name(),
                        runTime.timeName(),
                        mesh_
                    ),
                    mesh_,
sergio's avatar
sergio committed
664
                    dimensionedScalar("zero", dimensionSet(0, 3, -1, 0, 0), 0)
665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680
                )
            );
        }

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

sergio's avatar
sergio committed
681
            forAll(phases(), phasei)
682
            {
sergio's avatar
sergio committed
683
                alphaPhiSums[phasei] += phases()[phasei].alphaPhi();
684 685 686
            }
        }

sergio's avatar
sergio committed
687
        forAll(phases(), phasei)
688
        {
sergio's avatar
sergio committed
689
            phaseModel& phase = phases()[phasei];
690
            if (phase.stationary()) continue;
691

692
            phase.alphaPhiRef() = alphaPhiSums[phasei]/nAlphaSubCycles;
693 694 695 696 697 698 699
        }
    }
    else
    {
        solveAlphas();
    }

sergio's avatar
sergio committed
700
    forAll(phases(), phasei)
701
    {
sergio's avatar
sergio committed
702
        phaseModel& phase = phases()[phasei];
703 704 705 706
        if (phase.stationary()) continue;

        phase.alphaRhoPhiRef() =
            fvc::interpolate(phase.rho())*phase.alphaPhi();
707

708
        phase.clip(SMALL, 1 - SMALL);
709 710 711 712 713 714 715
    }

    calcAlphas();
}


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