twoPhaseSystem.C 8.99 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-2018 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 34
-------------------------------------------------------------------------------
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"
35
#include "UniformField.H"
36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63

#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),
64 65
    phase1_(phaseModels_[0]),
    phase2_(phaseModels_[1])
66 67
{
    phase2_.volScalarField::operator=(scalar(1) - phase1_);
68 69 70

    volScalarField& alpha1 = phase1_;
    mesh.setFluxRequired(alpha1.name());
71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123
}


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

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


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

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


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


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


void Foam::twoPhaseSystem::solve()
{
124
    const Time& runTime = mesh_.time();
125 126 127 128

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

129
    const dictionary& alphaControls = mesh_.solverDict(alpha1.name());
130

131 132
    label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
    label nAlphaCorr(alphaControls.get<label>("nAlphaCorr"));
133

134
    bool LTS = fv::localEulerDdt::enabled(mesh_);
135

136 137 138
    word alphaScheme("div(phi," + alpha1.name() + ')');
    word alpharScheme("div(phir," + alpha1.name() + ')');

139 140 141 142
    const surfaceScalarField& phi1 = phase1_.phi();
    const surfaceScalarField& phi2 = phase2_.phi();

    // Construct the dilatation rate source term
143
    tmp<volScalarField::Internal> tdgdt;
144

145
    if (phase1_.divU().valid() && phase2_.divU().valid())
146 147 148
    {
        tdgdt =
        (
149 150 151 152
            alpha2()
           *phase1_.divU()()()
          - alpha1()
           *phase2_.divU()()()
153 154
        );
    }
155
    else if (phase1_.divU().valid())
156 157 158
    {
        tdgdt =
        (
159 160
            alpha2()
           *phase1_.divU()()()
161 162
        );
    }
163
    else if (phase2_.divU().valid())
164 165 166
    {
        tdgdt =
        (
167 168
          - alpha1()
           *phase2_.divU()()()
169 170
        );
    }
171

172
    alpha1.correctBoundaryConditions();
173 174 175

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

176
    tmp<surfaceScalarField> alphaDbyA;
177
    if (DByAfs().found(phase1_.name()) && DByAfs().found(phase2_.name()))
178
    {
179 180 181 182
        surfaceScalarField DbyA
        (
            *DByAfs()[phase1_.name()] + *DByAfs()[phase2_.name()]
        );
183

184 185 186 187
        alphaDbyA =
            fvc::interpolate(max(alpha1, scalar(0)))
           *fvc::interpolate(max(alpha2, scalar(0)))
           *DbyA;
188

189
        phir += DbyA*fvc::snGrad(alpha1, "bounded")*mesh_.magSf();
190 191 192 193
    }

    for (int acorr=0; acorr<nAlphaCorr; acorr++)
    {
194
        volScalarField::Internal Sp
195 196 197 198 199
        (
            IOobject
            (
                "Sp",
                runTime.timeName(),
200
                mesh_
201
            ),
202
            mesh_,
203
            dimensionedScalar(dimless/dimTime)
204 205
        );

206
        volScalarField::Internal Su
207 208 209 210 211
        (
            IOobject
            (
                "Su",
                runTime.timeName(),
212
                mesh_
213 214 215
            ),
            // Divergence term is handled explicitly to be
            // consistent with the explicit transport solution
216
            fvc::div(phi_)*min(alpha1, scalar(1))
217 218
        );

219
        if (tdgdt.valid())
220
        {
221
            scalarField& dgdt = tdgdt.ref();
222 223

            forAll(dgdt, celli)
224
            {
225 226
                if (dgdt[celli] > 0.0)
                {
227 228
                    Sp[celli] -= dgdt[celli]/max(1 - alpha1[celli], 1e-4);
                    Su[celli] += dgdt[celli]/max(1 - alpha1[celli], 1e-4);
229 230 231 232 233
                }
                else if (dgdt[celli] < 0.0)
                {
                    Sp[celli] += dgdt[celli]/max(alpha1[celli], 1e-4);
                }
234 235 236
            }
        }

237
        surfaceScalarField alphaPhi1
238 239 240
        (
            fvc::flux
            (
241
                phi_,
242 243 244 245 246 247 248 249 250 251 252
                alpha1,
                alphaScheme
            )
          + fvc::flux
            (
               -fvc::flux(-phir, scalar(1) - alpha1, alpharScheme),
                alpha1,
                alpharScheme
            )
        );

253
        phase1_.correctInflowOutflow(alphaPhi1);
254 255 256

        if (nAlphaSubCycles > 1)
        {
257 258 259 260
            tmp<volScalarField> trSubDeltaT;

            if (LTS)
            {
261
                trSubDeltaT =
262
                    fv::localEulerDdt::localRSubDeltaT(mesh_, nAlphaSubCycles);
263 264
            }

265 266 267 268 269 270
            for
            (
                subCycle<volScalarField> alphaSubCycle(alpha1, nAlphaSubCycles);
                !(++alphaSubCycle).end();
            )
            {
271
                surfaceScalarField alphaPhi10(alphaPhi1);
272 273 274 275 276

                MULES::explicitSolve
                (
                    geometricOneField(),
                    alpha1,
277 278
                    phi_,
                    alphaPhi10,
279 280
                    (alphaSubCycle.index()*Sp)(),
                    (Su - (alphaSubCycle.index() - 1)*Sp*alpha1)(),
281 282
                    UniformField<scalar>(phase1_.alphaMax()),
                    zeroField()
283 284 285 286
                );

                if (alphaSubCycle.index() == 1)
                {
287
                    phase1_.alphaPhiRef() = alphaPhi10;
288 289 290
                }
                else
                {
291
                    phase1_.alphaPhiRef() += alphaPhi10;
292 293 294
                }
            }

295
            phase1_.alphaPhiRef() /= nAlphaSubCycles;
296 297 298 299 300 301 302
        }
        else
        {
            MULES::explicitSolve
            (
                geometricOneField(),
                alpha1,
303 304
                phi_,
                alphaPhi1,
305 306
                Sp,
                Su,
307 308
                UniformField<scalar>(phase1_.alphaMax()),
                zeroField()
309 310
            );

311
            phase1_.alphaPhiRef() = alphaPhi1;
312 313
        }

314
        if (alphaDbyA.valid())
315 316 317 318
        {
            fvScalarMatrix alpha1Eqn
            (
                fvm::ddt(alpha1) - fvc::ddt(alpha1)
319
              - fvm::laplacian(alphaDbyA(), alpha1, "bounded")
320 321 322 323 324
            );

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

325
            phase1_.alphaPhiRef() += alpha1Eqn.flux();
326 327
        }

328
        phase1_.alphaRhoPhiRef() =
329 330
            fvc::interpolate(phase1_.rho())*phase1_.alphaPhi();

331 332 333
        phase2_.alphaPhiRef() = phi_ - phase1_.alphaPhi();
        phase2_.correctInflowOutflow(phase2_.alphaPhiRef());
        phase2_.alphaRhoPhiRef() =
334 335 336
            fvc::interpolate(phase2_.rho())*phase2_.alphaPhi();

        Info<< alpha1.name() << " volume fraction = "
337
            << alpha1.weightedAverage(mesh_.V()).value()
338 339 340
            << "  Min(alpha1) = " << min(alpha1).value()
            << "  Max(alpha1) = " << max(alpha1).value()
            << endl;
341 342

        // Ensure the phase-fractions are bounded
343
        alpha1.clip(SMALL, 1 - SMALL);
344 345 346

        // Update the phase-fraction of the other phase
        alpha2 = scalar(1) - alpha1;
347 348 349 350 351
    }
}


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