twoPhaseSystem.C 10.5 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2013-2016 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 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59
     \\/     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 "twoPhaseSystem.H"
#include "dragModel.H"
#include "virtualMassModel.H"

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

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

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

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

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


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

Foam::twoPhaseSystem::twoPhaseSystem
(
    const fvMesh& mesh
)
:
    phaseSystem(mesh),
60 61
    phase1_(phaseModels_[0]),
    phase2_(phaseModels_[1])
62 63
{
    phase2_.volScalarField::operator=(scalar(1) - phase1_);
64 65 66

    volScalarField& alpha1 = phase1_;
    mesh.setFluxRequired(alpha1.name());
67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87
}


// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //

Foam::twoPhaseSystem::~twoPhaseSystem()
{}


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

Foam::tmp<Foam::volScalarField>
Foam::twoPhaseSystem::sigma() const
{
    return sigma
    (
        phasePairKey(phase1().name(), phase2().name())
    );
}


88 89 90 91 92 93
const Foam::dragModel& Foam::twoPhaseSystem::drag(const phaseModel& phase) const
{
    return lookupSubModel<dragModel>(phase, otherPhase(phase));
}


94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113
Foam::tmp<Foam::volScalarField>
Foam::twoPhaseSystem::Kd() const
{
    return Kd
    (
        phasePairKey(phase1().name(), phase2().name())
    );
}


Foam::tmp<Foam::surfaceScalarField>
Foam::twoPhaseSystem::Kdf() const
{
    return Kdf
    (
        phasePairKey(phase1().name(), phase2().name())
    );
}


114 115 116 117 118 119 120
const Foam::virtualMassModel&
Foam::twoPhaseSystem::virtualMass(const phaseModel& phase) const
{
    return lookupSubModel<virtualMassModel>(phase, otherPhase(phase));
}


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
Foam::tmp<Foam::volScalarField>
Foam::twoPhaseSystem::Vm() const
{
    return Vm
    (
        phasePairKey(phase1().name(), phase2().name())
    );
}


Foam::tmp<Foam::surfaceScalarField>
Foam::twoPhaseSystem::Vmf() const
{
    return Vmf
    (
        phasePairKey(phase1().name(), phase2().name())
    );
}


Foam::tmp<Foam::volVectorField>
Foam::twoPhaseSystem::F() const
{
    return F
    (
        phasePairKey(phase1().name(), phase2().name())
    );
}


Foam::tmp<Foam::surfaceScalarField>
Foam::twoPhaseSystem::Ff() const
{
    return Ff
    (
        phasePairKey(phase1().name(), phase2().name())
    );
}


Foam::tmp<Foam::volScalarField>
Foam::twoPhaseSystem::D() const
{
    return D
    (
        phasePairKey(phase1().name(), phase2().name())
    );
}


171 172 173 174 175 176
bool Foam::twoPhaseSystem::transfersMass() const
{
    return transfersMass(phase1());
}


177 178 179 180 181 182 183 184 185 186 187 188
Foam::tmp<Foam::volScalarField>
Foam::twoPhaseSystem::dmdt() const
{
    return dmdt
    (
        phasePairKey(phase1().name(), phase2().name())
    );
}


void Foam::twoPhaseSystem::solve()
{
189
    const Time& runTime = mesh_.time();
190 191 192 193

    volScalarField& alpha1 = phase1_;
    volScalarField& alpha2 = phase2_;

194
    const dictionary& alphaControls = mesh_.solverDict(alpha1.name());
195 196 197 198

    label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")));
    label nAlphaCorr(readLabel(alphaControls.lookup("nAlphaCorr")));

199
    bool LTS = fv::localEulerDdt::enabled(mesh_);
200

201 202 203
    word alphaScheme("div(phi," + alpha1.name() + ')');
    word alpharScheme("div(phir," + alpha1.name() + ')');

204 205 206 207 208
    const surfaceScalarField& phi = this->phi();
    const surfaceScalarField& phi1 = phase1_.phi();
    const surfaceScalarField& phi2 = phase2_.phi();

    // Construct the dilatation rate source term
209
    tmp<volScalarField::Internal> tdgdt;
210

211
    if (phase1_.divU().valid() && phase2_.divU().valid())
212 213 214
    {
        tdgdt =
        (
215 216 217 218
            alpha2()
           *phase1_.divU()()()
          - alpha1()
           *phase2_.divU()()()
219 220
        );
    }
221
    else if (phase1_.divU().valid())
222 223 224
    {
        tdgdt =
        (
225 226
            alpha2()
           *phase1_.divU()()()
227 228
        );
    }
229
    else if (phase2_.divU().valid())
230 231 232
    {
        tdgdt =
        (
233 234
          - alpha1()
           *phase2_.divU()()()
235 236
        );
    }
237

238 239
    alpha1.correctBoundaryConditions();
    surfaceScalarField alpha1f(fvc::interpolate(max(alpha1, scalar(0))));
240 241 242 243

    surfaceScalarField phic("phic", phi);
    surfaceScalarField phir("phir", phi1 - phi2);

244
    tmp<surfaceScalarField> alphaDbyA;
245

246
    if (notNull(phase1_.DbyA()) && notNull(phase2_.DbyA()))
247
    {
248
        surfaceScalarField DbyA(phase1_.DbyA() + phase2_.DbyA());
249

250 251 252 253
        alphaDbyA =
            fvc::interpolate(max(alpha1, scalar(0)))
           *fvc::interpolate(max(alpha2, scalar(0)))
           *DbyA;
254

255
        phir += DbyA*fvc::snGrad(alpha1, "bounded")*mesh_.magSf();
256 257 258 259
    }

    for (int acorr=0; acorr<nAlphaCorr; acorr++)
    {
260
        volScalarField::Internal Sp
261 262 263 264 265
        (
            IOobject
            (
                "Sp",
                runTime.timeName(),
266
                mesh_
267
            ),
268
            mesh_,
269
            dimensionedScalar("Sp", dimless/dimTime, 0.0)
270 271
        );

272
        volScalarField::Internal Su
273 274 275 276 277
        (
            IOobject
            (
                "Su",
                runTime.timeName(),
278
                mesh_
279 280 281 282 283 284
            ),
            // Divergence term is handled explicitly to be
            // consistent with the explicit transport solution
            fvc::div(phi)*min(alpha1, scalar(1))
        );

285
        if (tdgdt.valid())
286
        {
287
            scalarField& dgdt = tdgdt.ref();
288 289

            forAll(dgdt, celli)
290
            {
291 292 293 294 295 296 297 298 299
                if (dgdt[celli] > 0.0)
                {
                    Sp[celli] -= dgdt[celli]/max(1.0 - alpha1[celli], 1e-4);
                    Su[celli] += dgdt[celli]/max(1.0 - alpha1[celli], 1e-4);
                }
                else if (dgdt[celli] < 0.0)
                {
                    Sp[celli] += dgdt[celli]/max(alpha1[celli], 1e-4);
                }
300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318
            }
        }

        surfaceScalarField alphaPhic1
        (
            fvc::flux
            (
                phic,
                alpha1,
                alphaScheme
            )
          + fvc::flux
            (
               -fvc::flux(-phir, scalar(1) - alpha1, alpharScheme),
                alpha1,
                alpharScheme
            )
        );

319
        surfaceScalarField::Boundary& alphaPhic1Bf =
320 321
            alphaPhic1.boundaryFieldRef();

322
        // Ensure that the flux at inflow BCs is preserved
323
        forAll(alphaPhic1Bf, patchi)
324
        {
325
            fvsPatchScalarField& alphaPhic1p = alphaPhic1Bf[patchi];
326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343

            if (!alphaPhic1p.coupled())
            {
                const scalarField& phi1p = phi1.boundaryField()[patchi];
                const scalarField& alpha1p = alpha1.boundaryField()[patchi];

                forAll(alphaPhic1p, facei)
                {
                    if (phi1p[facei] < 0)
                    {
                        alphaPhic1p[facei] = alpha1p[facei]*phi1p[facei];
                    }
                }
            }
        }

        if (nAlphaSubCycles > 1)
        {
344 345 346 347
            tmp<volScalarField> trSubDeltaT;

            if (LTS)
            {
348
                trSubDeltaT =
349
                    fv::localEulerDdt::localRSubDeltaT(mesh_, nAlphaSubCycles);
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
            for
            (
                subCycle<volScalarField> alphaSubCycle(alpha1, nAlphaSubCycles);
                !(++alphaSubCycle).end();
            )
            {
                surfaceScalarField alphaPhic10(alphaPhic1);

                MULES::explicitSolve
                (
                    geometricOneField(),
                    alpha1,
                    phi,
                    alphaPhic10,
                    (alphaSubCycle.index()*Sp)(),
                    (Su - (alphaSubCycle.index() - 1)*Sp*alpha1)(),
                    phase1_.alphaMax(),
                    0
                );

                if (alphaSubCycle.index() == 1)
                {
                    phase1_.alphaPhi() = alphaPhic10;
                }
                else
                {
                    phase1_.alphaPhi() += alphaPhic10;
                }
            }

            phase1_.alphaPhi() /= nAlphaSubCycles;
        }
        else
        {
            MULES::explicitSolve
            (
                geometricOneField(),
                alpha1,
                phi,
                alphaPhic1,
                Sp,
                Su,
                phase1_.alphaMax(),
                0
            );

            phase1_.alphaPhi() = alphaPhic1;
        }

401
        if (alphaDbyA.valid())
402 403 404 405
        {
            fvScalarMatrix alpha1Eqn
            (
                fvm::ddt(alpha1) - fvc::ddt(alpha1)
406
              - fvm::laplacian(alphaDbyA, alpha1, "bounded")
407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423
            );

            alpha1Eqn.relax();
            alpha1Eqn.solve();

            phase1_.alphaPhi() += alpha1Eqn.flux();
        }

        phase1_.alphaRhoPhi() =
            fvc::interpolate(phase1_.rho())*phase1_.alphaPhi();

        phase2_.alphaPhi() = phi - phase1_.alphaPhi();
        alpha2 = scalar(1) - alpha1;
        phase2_.alphaRhoPhi() =
            fvc::interpolate(phase2_.rho())*phase2_.alphaPhi();

        Info<< alpha1.name() << " volume fraction = "
424
            << alpha1.weightedAverage(mesh_.V()).value()
425 426 427 428 429 430 431 432
            << "  Min(alpha1) = " << min(alpha1).value()
            << "  Max(alpha1) = " << max(alpha1).value()
            << endl;
    }
}


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