pEqn.H 9.7 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
surfaceScalarField alphaf1("alphaf1", fvc::interpolate(alpha1));
surfaceScalarField alphaf2("alphaf2", scalar(1) - alphaf1);

surfaceScalarField alphaRhof10
(
    "alphaRhof10",
    fvc::interpolate
    (
        max(alpha1.oldTime(), phase1.residualAlpha())
       *rho1.oldTime()
    )
);

surfaceScalarField alphaRhof20
(
    "alphaRhof20",
    fvc::interpolate
    (
        max(alpha2.oldTime(), phase2.residualAlpha())
       *rho2.oldTime()
    )
);

// Drag coefficient
surfaceScalarField Kdf("Kdf", fluid.Kdf());

// Virtual-mass coefficient
surfaceScalarField Vmf("Vmf", fluid.Vmf());

surfaceScalarField rAUf1
(
    IOobject::groupName("rAUf", phase1.name()),
    1.0
   /(
35
        byDt(alphaRhof10 + Vmf)
36
37
      + fvc::interpolate(U1Eqn.A())
      + Kdf
38
    )
39
40
41
42
43
44
45
);

surfaceScalarField rAUf2
(
    IOobject::groupName("rAUf", phase2.name()),
    1.0
   /(
46
        byDt(alphaRhof20 + Vmf)
47
48
      + fvc::interpolate(U2Eqn.A())
      + Kdf
49
    )
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
);


// Turbulent dispersion, particle-pressure, lift and wall-lubrication forces
tmp<surfaceScalarField> Ff1;
tmp<surfaceScalarField> Ff2;
{
    // Turbulent-dispersion diffusivity
    volScalarField D(fluid.D());

    // Phase-1 turbulent dispersion and particle-pressure diffusivity
    surfaceScalarField Df1
    (
        fvc::interpolate(D + phase1.turbulence().pPrime())
    );

    // Phase-2 turbulent dispersion and particle-pressure diffusivity
    surfaceScalarField Df2
    (
        fvc::interpolate(D + phase2.turbulence().pPrime())
    );

72
    // Cache the phase diffusivities for implicit treatment in the
73
74
75
    // phase-fraction equation
    if (implicitPhasePressure)
    {
76
77
        phase1.DbyA(rAUf1*Df1);
        phase2.DbyA(rAUf2*Df2);
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
    }

    // Lift and wall-lubrication forces
    surfaceScalarField Ff(fluid.Ff());

    // Phase-fraction face-gradient
    surfaceScalarField snGradAlpha1(fvc::snGrad(alpha1)*mesh.magSf());

    // Phase-1 dispersion, lift and wall-lubrication force
    Ff1 = Df1*snGradAlpha1 + Ff;

    // Phase-2 dispersion, lift and wall-lubrication force
    Ff2 = -Df2*snGradAlpha1 - Ff;
}


while (pimple.correct())
{
96
97
98
99
    volScalarField rho("rho", fluid.rho());

    // Correct p_rgh for consistency with p and the updated densities
    p_rgh = p - rho*gh;
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
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

    surfaceScalarField rhof1(fvc::interpolate(rho1));
    surfaceScalarField rhof2(fvc::interpolate(rho2));

    // Correct fixed-flux BCs to be consistent with the velocity BCs
    MRF.correctBoundaryFlux(U1, phi1);
    MRF.correctBoundaryFlux(U2, phi2);

    surfaceScalarField alpharAUf1
    (
        IOobject::groupName("alpharAUf", phase1.name()),
        max(alphaf1, phase1.residualAlpha())*rAUf1
    );

    surfaceScalarField alpharAUf2
    (
        IOobject::groupName("alpharAUf", phase2.name()),
        max(alphaf2, phase2.residualAlpha())*rAUf2
    );

    surfaceScalarField ghSnGradRho
    (
        "ghSnGradRho",
        ghf*fvc::snGrad(rho)*mesh.magSf()
    );

    // Phase-1 buoyancy flux
    surfaceScalarField phig1
    (
        IOobject::groupName("phig", phase1.name()),
        alpharAUf1
       *(
            ghSnGradRho
          - alphaf2*(rhof1 - rhof2)*(g & mesh.Sf())
        )
    );

    // Phase-2 buoyancy flux
    surfaceScalarField phig2
    (
        IOobject::groupName("phig", phase2.name()),
        alpharAUf2
       *(
            ghSnGradRho
          - alphaf1*(rhof2 - rhof1)*(g & mesh.Sf())
        )
    );


    // Phase-1 predicted flux
    surfaceScalarField phiHbyA1
    (
        IOobject::groupName("phiHbyA", phase1.name()),
        phi1
    );

    phiHbyA1 =
        rAUf1
       *(
            (alphaRhof10 + Vmf)
160
           *byDt(MRF.absolute(phi1.oldTime()))
161
          + fvc::flux(U1Eqn.H())
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
          + Vmf*ddtPhi2
          + Kdf*MRF.absolute(phi2)
          - Ff1()
        );

    // Phase-2 predicted flux
    surfaceScalarField phiHbyA2
    (
        IOobject::groupName("phiHbyA", phase2.name()),
        phi2
    );

    phiHbyA2 =
        rAUf2
       *(
            (alphaRhof20 + Vmf)
178
           *byDt(MRF.absolute(phi2.oldTime()))
179
          + fvc::flux(U2Eqn.H())
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
          + Vmf*ddtPhi1
          + Kdf*MRF.absolute(phi1)
          - Ff2()
       );


    surfaceScalarField phiHbyA
    (
        "phiHbyA",
        alphaf1*(phiHbyA1 - phig1) + alphaf2*(phiHbyA2 - phig2)
    );
    MRF.makeRelative(phiHbyA);

    phiHbyA1 -= phig1;
    phiHbyA2 -= phig2;

    surfaceScalarField rAUf
    (
        "rAUf",
        mag(alphaf1*alpharAUf1 + alphaf2*alpharAUf2)
    );

    // Update the fixedFluxPressure BCs to ensure flux consistency
    setSnGrad<fixedFluxPressureFvPatchScalarField>
    (
205
        p_rgh.boundaryFieldRef(),
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
        (
            phiHbyA.boundaryField()
          - (
                alphaf1.boundaryField()*phi1.boundaryField()
              + alphaf2.boundaryField()*phi2.boundaryField()
            )
        )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
    );

    tmp<fvScalarMatrix> pEqnComp1;
    tmp<fvScalarMatrix> pEqnComp2;

    if (pimple.transonic())
    {
        surfaceScalarField phid1
        (
            IOobject::groupName("phid", phase1.name()),
            fvc::interpolate(psi1)*phi1
        );
        surfaceScalarField phid2
        (
            IOobject::groupName("phid", phase2.name()),
            fvc::interpolate(psi2)*phi2
        );

231
232
233
        if (phase1.compressible())
        {
            pEqnComp1 =
234
            (
235
                phase1.continuityError()
236
237
              - fvc::Sp(fvc::ddt(alpha1) + fvc::div(alphaPhi1), rho1)
            )/rho1
238
          + correction
239
            (
240
241
242
243
244
                (alpha1/rho1)*
                (
                    psi1*fvm::ddt(p_rgh)
                  + fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
                )
245
            );
246
247
            deleteDemandDrivenData(pEqnComp1.ref().faceFluxCorrectionPtr());
            pEqnComp1.ref().relax();
248
        }
249

250
251
252
        if (phase2.compressible())
        {
            pEqnComp2 =
253
            (
254
                phase2.continuityError()
255
256
              - fvc::Sp(fvc::ddt(alpha2) + fvc::div(alphaPhi2), rho2)
            )/rho2
257
          + correction
258
            (
259
260
261
262
263
                (alpha2/rho2)*
                (
                    psi2*fvm::ddt(p_rgh)
                  + fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
                )
264
            );
265
266
            deleteDemandDrivenData(pEqnComp2.ref().faceFluxCorrectionPtr());
            pEqnComp2.ref().relax();
267
        }
268
269
270
    }
    else
    {
271
272
273
        if (phase1.compressible())
        {
            pEqnComp1 =
274
            (
275
                phase1.continuityError()
276
277
278
              - fvc::Sp(fvc::ddt(alpha1) + fvc::div(alphaPhi1), rho1)
            )/rho1
          + (alpha1*psi1/rho1)*correction(fvm::ddt(p_rgh));
279
        }
280

281
282
283
        if (phase2.compressible())
        {
            pEqnComp2 =
284
            (
285
                phase2.continuityError()
286
287
288
              - fvc::Sp(fvc::ddt(alpha2) + fvc::div(alphaPhi2), rho2)
            )/rho2
          + (alpha2*psi2/rho2)*correction(fvm::ddt(p_rgh));
289
        }
290
291
    }

292
293
294
295
    if (fluid.transfersMass())
    {
        if (pEqnComp1.valid())
        {
296
            pEqnComp1.ref() -= fluid.dmdt()/rho1;
297
298
299
300
301
302
303
304
        }
        else
        {
            pEqnComp1 = fvm::Su(-fluid.dmdt()/rho1, p_rgh);
        }

        if (pEqnComp2.valid())
        {
305
            pEqnComp2.ref() += fluid.dmdt()/rho2;
306
307
308
309
310
311
312
        }
        else
        {
            pEqnComp2 = fvm::Su(fluid.dmdt()/rho2, p_rgh);
        }
    }

313
314
315
316
317
318
319
320
321
322
323
    // Cache p prior to solve for density update
    volScalarField p_rgh_0("p_rgh_0", p_rgh);

    while (pimple.correctNonOrthogonal())
    {
        fvScalarMatrix pEqnIncomp
        (
            fvc::div(phiHbyA)
          - fvm::laplacian(rAUf, p_rgh)
        );

324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
        {
            fvScalarMatrix pEqn(pEqnIncomp);

            if (pEqnComp1.valid())
            {
                pEqn += pEqnComp1();
            }

            if (pEqnComp2.valid())
            {
                pEqn += pEqnComp2();
            }

            solve
            (
                pEqn,
                mesh.solver(p_rgh.select(pimple.finalInnerIter()))
            );
        }
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

        if (pimple.finalNonOrthogonalIter())
        {
            surfaceScalarField mSfGradp("mSfGradp", pEqnIncomp.flux()/rAUf);

            phi = phiHbyA + pEqnIncomp.flux();

            surfaceScalarField phi1s
            (
                phiHbyA1
              + alpharAUf1*mSfGradp
              - rAUf1*Kdf*MRF.absolute(phi2)
            );

            surfaceScalarField phi2s
            (
                phiHbyA2
              + alpharAUf2*mSfGradp
              - rAUf2*Kdf*MRF.absolute(phi1)
            );

            surfaceScalarField phir
            (
                ((phi2s + rAUf2*Kdf*phi1s) - (phi1s + rAUf1*Kdf*phi2s))
               /(1.0 - rAUf1*rAUf2*sqr(Kdf))
            );

            phi1 = phi - alphaf2*phir;
            phi2 = phi + alphaf1*phir;

            U1 = fvc::reconstruct(MRF.absolute(phi1));
            U1.correctBoundaryConditions();
            fvOptions.correct(U1);

            U2 = fvc::reconstruct(MRF.absolute(phi2));
            U2.correctBoundaryConditions();
            fvOptions.correct(U2);

381
            // Set the phase dilatation rates
382
            if (pEqnComp1.valid())
383
            {
384
                phase1.divU(-pEqnComp1 & p_rgh);
385
            }
386
            if (pEqnComp2.valid())
387
            {
388
                phase2.divU(-pEqnComp2 & p_rgh);
389
            }
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
        }
    }

    // Update and limit the static pressure
    p = max(p_rgh + rho*gh, pMin);

    // Limit p_rgh
    p_rgh = p - rho*gh;

    // Update densities from change in p_rgh
    rho1 += psi1*(p_rgh - p_rgh_0);
    rho2 += psi2*(p_rgh - p_rgh_0);

    // Correct p_rgh for consistency with p and the updated densities
    rho = fluid.rho();
    p_rgh = p - rho*gh;
    p_rgh.correctBoundaryConditions();
}