twoPhaseSystem.C 10.6 KB
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
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2013-2015 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/>.

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

#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_.first()),
    phase2_(phaseModels_.last())
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 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193
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())
    );
}


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


void Foam::twoPhaseSystem::solve()
{
    const fvMesh& mesh = this->mesh();
    const Time& runTime = mesh.time();

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

    const dictionary& alphaControls = mesh.solverDict(alpha1.name());

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

194
    bool LTS = fv::localEulerDdt::enabled(mesh);
195

196 197 198
    word alphaScheme("div(phi," + alpha1.name() + ')');
    word alpharScheme("div(phir," + alpha1.name() + ')');

199 200 201 202 203 204 205 206 207 208 209
    const surfaceScalarField& phi = this->phi();
    const surfaceScalarField& phi1 = phase1_.phi();
    const surfaceScalarField& phi2 = phase2_.phi();

    // Construct the dilatation rate source term
    tmp<volScalarField::DimensionedInternalField> tdgdt;

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

233 234
    alpha1.correctBoundaryConditions();
    surfaceScalarField alpha1f(fvc::interpolate(max(alpha1, scalar(0))));
235 236 237 238

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

239
    tmp<surfaceScalarField> alphaDbyA;
240

241
    if (notNull(phase1_.DbyA()) && notNull(phase2_.DbyA()))
242
    {
243
        surfaceScalarField DbyA(phase1_.DbyA() + phase2_.DbyA());
244

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

250
        phir += DbyA*fvc::snGrad(alpha1, "bounded")*mesh_.magSf();
251 252 253 254 255 256 257 258 259 260 261 262 263
    }

    for (int acorr=0; acorr<nAlphaCorr; acorr++)
    {
        volScalarField::DimensionedInternalField Sp
        (
            IOobject
            (
                "Sp",
                runTime.timeName(),
                mesh
            ),
            mesh,
264
            dimensionedScalar("Sp", dimless/dimTime, 0.0)
265 266 267 268 269 270 271 272 273 274 275 276 277 278 279
        );

        volScalarField::DimensionedInternalField Su
        (
            IOobject
            (
                "Su",
                runTime.timeName(),
                mesh
            ),
            // Divergence term is handled explicitly to be
            // consistent with the explicit transport solution
            fvc::div(phi)*min(alpha1, scalar(1))
        );

280
        if (tdgdt.valid())
281
        {
282 283 284
            scalarField& dgdt = tdgdt();

            forAll(dgdt, celli)
285
            {
286 287 288 289 290 291 292 293 294
                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);
                }
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
            }
        }

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

        // Ensure that the flux at inflow BCs is preserved
        forAll(alphaPhic1.boundaryField(), patchi)
        {
            fvsPatchScalarField& alphaPhic1p =
                alphaPhic1.boundaryField()[patchi];

            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)
        {
337 338 339 340
            tmp<volScalarField> trSubDeltaT;

            if (LTS)
            {
341 342
                trSubDeltaT =
                    fv::localEulerDdt::localRSubDeltaT(mesh, nAlphaSubCycles);
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
            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;
        }

394
        if (alphaDbyA.valid())
395 396 397 398
        {
            fvScalarMatrix alpha1Eqn
            (
                fvm::ddt(alpha1) - fvc::ddt(alpha1)
399
              - fvm::laplacian(alphaDbyA, alpha1, "bounded")
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
            );

            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 = "
            << alpha1.weightedAverage(mesh.V()).value()
            << "  Min(alpha1) = " << min(alpha1).value()
            << "  Max(alpha1) = " << max(alpha1).value()
            << endl;
    }
}


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