twoPhaseSystem.C 9.84 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 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
    label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
    label nAlphaCorr(alphaControls.get<label>("nAlphaCorr"));
198

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
    const surfaceScalarField& phi1 = phase1_.phi();
    const surfaceScalarField& phi2 = phase2_.phi();

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

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

237
    alpha1.correctBoundaryConditions();
238 239 240

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

241
    tmp<surfaceScalarField> alphaDbyA;
242

243
    if (notNull(phase1_.DbyA()) && notNull(phase2_.DbyA()))
244
    {
245
        surfaceScalarField DbyA(phase1_.DbyA() + phase2_.DbyA());
246

247 248 249 250
        alphaDbyA =
            fvc::interpolate(max(alpha1, scalar(0)))
           *fvc::interpolate(max(alpha2, scalar(0)))
           *DbyA;
251

252
        phir += DbyA*fvc::snGrad(alpha1, "bounded")*mesh_.magSf();
253 254 255 256
    }

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

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

282
        if (tdgdt.valid())
283
        {
284
            scalarField& dgdt = tdgdt.ref();
285 286

            forAll(dgdt, celli)
287
            {
288 289 290 291 292 293 294 295 296
                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);
                }
297 298 299
            }
        }

300
        surfaceScalarField alphaPhi1
301 302 303
        (
            fvc::flux
            (
304
                phi_,
305 306 307 308 309 310 311 312 313 314 315
                alpha1,
                alphaScheme
            )
          + fvc::flux
            (
               -fvc::flux(-phir, scalar(1) - alpha1, alpharScheme),
                alpha1,
                alpharScheme
            )
        );

316
        phase1_.correctInflowOutflow(alphaPhi1);
317 318 319

        if (nAlphaSubCycles > 1)
        {
320 321 322 323
            tmp<volScalarField> trSubDeltaT;

            if (LTS)
            {
324
                trSubDeltaT =
325
                    fv::localEulerDdt::localRSubDeltaT(mesh_, nAlphaSubCycles);
326 327
            }

328 329 330 331 332 333
            for
            (
                subCycle<volScalarField> alphaSubCycle(alpha1, nAlphaSubCycles);
                !(++alphaSubCycle).end();
            )
            {
334
                surfaceScalarField alphaPhi10(alphaPhi1);
335 336 337 338 339

                MULES::explicitSolve
                (
                    geometricOneField(),
                    alpha1,
340 341
                    phi_,
                    alphaPhi10,
342 343 344 345 346 347 348 349
                    (alphaSubCycle.index()*Sp)(),
                    (Su - (alphaSubCycle.index() - 1)*Sp*alpha1)(),
                    phase1_.alphaMax(),
                    0
                );

                if (alphaSubCycle.index() == 1)
                {
350
                    phase1_.alphaPhi() = alphaPhi10;
351 352 353
                }
                else
                {
354
                    phase1_.alphaPhi() += alphaPhi10;
355 356 357 358 359 360 361 362 363 364 365
                }
            }

            phase1_.alphaPhi() /= nAlphaSubCycles;
        }
        else
        {
            MULES::explicitSolve
            (
                geometricOneField(),
                alpha1,
366 367
                phi_,
                alphaPhi1,
368 369 370 371 372 373
                Sp,
                Su,
                phase1_.alphaMax(),
                0
            );

374
            phase1_.alphaPhi() = alphaPhi1;
375 376
        }

377
        if (alphaDbyA.valid())
378 379 380 381
        {
            fvScalarMatrix alpha1Eqn
            (
                fvm::ddt(alpha1) - fvc::ddt(alpha1)
382
              - fvm::laplacian(alphaDbyA(), alpha1, "bounded")
383 384 385 386 387 388 389 390 391 392 393
            );

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

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

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

394 395
        phase2_.alphaPhi() = phi_ - phase1_.alphaPhi();
        phase2_.correctInflowOutflow(phase2_.alphaPhi());
396 397 398 399
        phase2_.alphaRhoPhi() =
            fvc::interpolate(phase2_.rho())*phase2_.alphaPhi();

        Info<< alpha1.name() << " volume fraction = "
400
            << alpha1.weightedAverage(mesh_.V()).value()
401 402 403
            << "  Min(alpha1) = " << min(alpha1).value()
            << "  Max(alpha1) = " << max(alpha1).value()
            << endl;
404 405

        // Ensure the phase-fractions are bounded
406
        alpha1.clip(0, 1);
407 408 409

        // Update the phase-fraction of the other phase
        alpha2 = scalar(1) - alpha1;
410 411 412 413 414
    }
}


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