pEqn.H
Go to the documentation of this file.
2 surfaceScalarField alphaf2("alphaf2", scalar(1) - alphaf1);
3 
5 (
6  IOobject::groupName("rAU", phase1.name()),
7  1.0
8  /(
9  U1Eqn.A()
10  + max(phase1.residualAlpha() - alpha1, scalar(0))
11  *rho1/runTime.deltaT()
12  )
13 );
15 (
16  IOobject::groupName("rAU", phase2.name()),
17  1.0
18  /(
19  U2Eqn.A()
20  + max(phase2.residualAlpha() - alpha2, scalar(0))
21  *rho2/runTime.deltaT()
22  )
23 );
24 
26 (
27  fvc::interpolate(max(alpha1, phase1.residualAlpha())*rAU1)
28 );
30 (
31  fvc::interpolate(max(alpha2, phase2.residualAlpha())*rAU2)
32 );
33 
34 volScalarField Kd(fluid.Kd());
35 
36 // Turbulent diffusion, particle-pressure, lift and wall-lubrication fluxes
37 tmp<surfaceScalarField> phiF1;
38 tmp<surfaceScalarField> phiF2;
39 {
40  // Turbulent-dispersion diffusivity
41  volScalarField D(fluid.D());
42 
43  // Phase-1 turbulent dispersion and particle-pressure flux
44  tmp<surfaceScalarField> DbyA1
45  (
47  (
48  rAU1*(D + phase1.turbulence().pPrime())
49  )
50  );
51 
52  // Phase-2 turbulent dispersion and particle-pressure flux
53  tmp<surfaceScalarField> DbyA2
54  (
56  (
57  rAU2*(D + phase2.turbulence().pPrime())
58  )
59  );
60 
61  // Lift and wall-lubrication forces
62  volVectorField F(fluid.F());
63 
64  // Phase-fraction face-gradient
66 
67  // Phase-1 dispersion, lift and wall-lubrication flux
68  phiF1 = DbyA1()*snGradAlpha1 + fvc::flux(rAU1*F);
69 
70  // Phase-2 dispersion, lift and wall-lubrication flux
71  phiF2 = - DbyA2()*snGradAlpha1 - fvc::flux(rAU2*F);
72 
73  // Cache the phase diffusivities for implicit treatment in the
74  // phase-fraction equation
75  if (implicitPhasePressure)
76  {
77  phase1.DbyA(DbyA1);
78  phase2.DbyA(DbyA2);
79  }
80 }
81 
82 
83 // --- Pressure corrector loop
84 while (pimple.correct())
85 {
86  // Update continuity errors due to temperature changes
87  fluid.correct();
88 
89  volScalarField rho("rho", fluid.rho());
90 
91  // Correct p_rgh for consistency with p and the updated densities
92  p_rgh = p - rho*gh;
93 
94  // Correct fixed-flux BCs to be consistent with the velocity BCs
95  MRF.correctBoundaryFlux(U1, phi1);
96  MRF.correctBoundaryFlux(U2, phi2);
97 
98  volVectorField HbyA1
99  (
100  IOobject::groupName("HbyA", phase1.name()),
101  U1
102  );
103  HbyA1 =
104  rAU1
105  *(
106  U1Eqn.H()
107  + max(phase1.residualAlpha() - alpha1, scalar(0))
108  *rho1*U1.oldTime()/runTime.deltaT()
109  );
110 
111  volVectorField HbyA2
112  (
113  IOobject::groupName("HbyA", phase2.name()),
114  U2
115  );
116  HbyA2 =
117  rAU2
118  *(
119  U2Eqn.H()
120  + max(phase2.residualAlpha() - alpha2, scalar(0))
121  *rho2*U2.oldTime()/runTime.deltaT()
122  );
123 
125  (
126  "ghSnGradRho",
127  ghf*fvc::snGrad(rho)*mesh.magSf()
128  );
129 
130  surfaceScalarField phig1
131  (
132  alpharAUf1
133  *(
135  - alphaf2*fvc::interpolate(rho1 - rho2)*(g & mesh.Sf())
136  )
137  );
138 
139  surfaceScalarField phig2
140  (
141  alpharAUf2
142  *(
144  - alphaf1*fvc::interpolate(rho2 - rho1)*(g & mesh.Sf())
145  )
146  );
147 
148 
149  // ddtPhiCorr filter -- only apply in pure(ish) phases
151  surfaceScalarField phiCorrCoeff1(pos(alphaf1Bar - 0.99));
152  surfaceScalarField phiCorrCoeff2(pos(0.01 - alphaf1Bar));
153 
154  {
155  surfaceScalarField::Boundary& phiCorrCoeff1Bf =
156  phiCorrCoeff1.boundaryFieldRef();
157 
158  surfaceScalarField::Boundary& phiCorrCoeff2Bf =
159  phiCorrCoeff2.boundaryFieldRef();
160 
161  forAll(mesh.boundary(), patchi)
162  {
163  // Set ddtPhiCorr to 0 on non-coupled boundaries
164  if
165  (
166  !mesh.boundary()[patchi].coupled()
167  || isA<cyclicAMIFvPatch>(mesh.boundary()[patchi])
168  )
169  {
170  phiCorrCoeff1Bf[patchi] = 0;
171  phiCorrCoeff2Bf[patchi] = 0;
172  }
173  }
174  }
175 
176  // Phase-1 predicted flux
177  surfaceScalarField phiHbyA1
178  (
179  IOobject::groupName("phiHbyA", phase1.name()),
180  fvc::flux(HbyA1)
181  + phiCorrCoeff1*fvc::interpolate(alpha1.oldTime()*rho1.oldTime()*rAU1)
182  *(
183  MRF.absolute(phi1.oldTime())
184  - fvc::flux(U1.oldTime())
185  )/runTime.deltaT()
186  - phiF1()
187  - phig1
188  );
189 
190  // Phase-2 predicted flux
191  surfaceScalarField phiHbyA2
192  (
193  IOobject::groupName("phiHbyA", phase2.name()),
194  fvc::flux(HbyA2)
195  + phiCorrCoeff2*fvc::interpolate(alpha2.oldTime()*rho2.oldTime()*rAU2)
196  *(
197  MRF.absolute(phi2.oldTime())
198  - fvc::flux(U2.oldTime())
199  )/runTime.deltaT()
200  - phiF2()
201  - phig2
202  );
203 
204  // Face-drag coefficients
207 
208  // Construct the mean predicted flux
209  // including explicit drag contributions based on absolute fluxes
211  (
212  "phiHbyA",
213  alphaf1*(phiHbyA1 + rAUKd1*MRF.absolute(phi2))
214  + alphaf2*(phiHbyA2 + rAUKd2*MRF.absolute(phi1))
215  );
216  MRF.makeRelative(phiHbyA);
217 
218  // Construct pressure "diffusivity"
220  (
221  "rAUf",
223  );
224 
225  // Update the fixedFluxPressure BCs to ensure flux consistency
227  (
228  p_rgh.boundaryFieldRef(),
229  (
230  phiHbyA.boundaryField()
231  - (
232  alphaf1.boundaryField()*phi1.boundaryField()
233  + alphaf2.boundaryField()*phi2.boundaryField()
234  )
235  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
236  );
237 
238  tmp<fvScalarMatrix> pEqnComp1;
239  tmp<fvScalarMatrix> pEqnComp2;
240 
241  // Construct the compressibility parts of the pressure equation
242  if (pimple.transonic())
243  {
244  if (phase1.compressible())
245  {
246  surfaceScalarField phid1
247  (
248  IOobject::groupName("phid", phase1.name()),
250  );
251 
252  pEqnComp1 =
253  (
254  phase1.continuityError()
256  )/rho1
257  + correction
258  (
259  (alpha1/rho1)*
260  (
261  psi1*fvm::ddt(p_rgh)
262  + fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
263  )
264  );
265 
266  deleteDemandDrivenData(pEqnComp1.ref().faceFluxCorrectionPtr());
267  pEqnComp1.ref().relax();
268  }
269 
270  if (phase2.compressible())
271  {
272  surfaceScalarField phid2
273  (
274  IOobject::groupName("phid", phase2.name()),
276  );
277 
278  pEqnComp2 =
279  (
280  phase2.continuityError()
282  )/rho2
283  + correction
284  (
285  (alpha2/rho2)*
286  (
287  psi2*fvm::ddt(p_rgh)
288  + fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
289  )
290  );
291  deleteDemandDrivenData(pEqnComp2.ref().faceFluxCorrectionPtr());
292  pEqnComp2.ref().relax();
293  }
294  }
295  else
296  {
297  if (phase1.compressible())
298  {
299  pEqnComp1 =
300  (
301  phase1.continuityError()
303  )/rho1
304  + (alpha1*psi1/rho1)*correction(fvm::ddt(p_rgh));
305  }
306 
307  if (phase2.compressible())
308  {
309  pEqnComp2 =
310  (
311  phase2.continuityError()
313  )/rho2
314  + (alpha2*psi2/rho2)*correction(fvm::ddt(p_rgh));
315  }
316  }
317 
318  if (fluid.transfersMass())
319  {
320  if (pEqnComp1.valid())
321  {
322  pEqnComp1.ref() -= fluid.dmdt()/rho1;
323  }
324  else
325  {
326  pEqnComp1 = fvm::Su(-fluid.dmdt()/rho1, p_rgh);
327  }
328 
329  if (pEqnComp2.valid())
330  {
331  pEqnComp2.ref() += fluid.dmdt()/rho2;
332  }
333  else
334  {
335  pEqnComp2 = fvm::Su(fluid.dmdt()/rho2, p_rgh);
336  }
337  }
338 
339  // Cache p prior to solve for density update
340  volScalarField p_rgh_0(p_rgh);
341 
342  // Iterate over the pressure equation to correct for non-orthogonality
343  while (pimple.correctNonOrthogonal())
344  {
345  // Construct the transport part of the pressure equation
346  fvScalarMatrix pEqnIncomp
347  (
348  fvc::div(phiHbyA)
349  - fvm::laplacian(rAUf, p_rgh)
350  );
351 
352  {
353  fvScalarMatrix pEqn(pEqnIncomp);
354 
355  if (pEqnComp1.valid())
356  {
357  pEqn += pEqnComp1();
358  }
359 
360  if (pEqnComp2.valid())
361  {
362  pEqn += pEqnComp2();
363  }
364 
365  solve
366  (
367  pEqn,
368  mesh.solver(p_rgh.select(pimple.finalInnerIter()))
369  );
370  }
371 
372  // Correct fluxes and velocities on last non-orthogonal iteration
373  if (pimple.finalNonOrthogonalIter())
374  {
375  phi = phiHbyA + pEqnIncomp.flux();
376 
377  surfaceScalarField mSfGradp("mSfGradp", pEqnIncomp.flux()/rAUf);
378 
379  // Partial-elimination phase-flux corrector
380  {
381  surfaceScalarField phi1s
382  (
383  phiHbyA1 + alpharAUf1*mSfGradp
384  );
385 
386  surfaceScalarField phi2s
387  (
388  phiHbyA2 + alpharAUf2*mSfGradp
389  );
390 
392  (
393  ((phi1s + rAUKd1*phi2s) - (phi2s + rAUKd2*phi1s))
394  /(1 - rAUKd1*rAUKd2)
395  );
396 
397  phi1 = phi + alphaf2*phir;
398  phi2 = phi - alphaf1*phir;
399  }
400 
401  // Set the phase dilatation rates
402  if (phase1.compressible())
403  {
404  phase1.divU(-pEqnComp1 & p_rgh);
405  }
406  if (phase2.compressible())
407  {
408  phase2.divU(-pEqnComp2 & p_rgh);
409  }
410 
411  // Optionally relax pressure for velocity correction
412  p_rgh.relax();
413 
414  mSfGradp = pEqnIncomp.flux()/rAUf;
415 
416  // Partial-elimination phase-velocity corrector
417  {
418  volVectorField Us1
419  (
420  HbyA1
421  + fvc::reconstruct(alpharAUf1*mSfGradp - phiF1() - phig1)
422  );
423 
424  volVectorField Us2
425  (
426  HbyA2
427  + fvc::reconstruct(alpharAUf2*mSfGradp - phiF2() - phig2)
428  );
429 
430  volScalarField D1(rAU1*Kd);
431  volScalarField D2(rAU2*Kd);
432 
433  volVectorField U(alpha1*(Us1 + D1*U2) + alpha2*(Us2 + D2*U1));
434  volVectorField Ur(((1 - D2)*Us1 - (1 - D1)*Us2)/(1 - D1*D2));
435 
436  U1 = U + alpha2*Ur;
437  U1.correctBoundaryConditions();
438  fvOptions.correct(U1);
439 
440  U2 = U - alpha1*Ur;
441  U2.correctBoundaryConditions();
442  fvOptions.correct(U2);
443  }
444  }
445  }
446 
447  // Update and limit the static pressure
448  p = max(p_rgh + rho*gh, pMin);
449 
450  // Limit p_rgh
451  p_rgh = p - rho*gh;
452 
453  // Update densities from change in p_rgh
454  rho1 += psi1*(p_rgh - p_rgh_0);
455  rho2 += psi2*(p_rgh - p_rgh_0);
456 
457  // Correct p_rgh for consistency with p and the updated densities
458  rho = fluid.rho();
459  p_rgh = p - rho*gh;
460  p_rgh.correctBoundaryConditions();
461 }
volScalarField rAU1(IOobject::groupName("rAU", phase1.name()), 1.0/(U1Eqn.A()+max(phase1.residualAlpha()-alpha1, scalar(0))*rho1/runTime.deltaT()))
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
psi2
Definition: TEqns.H:35
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
tmp< surfaceScalarField > phiF2
Definition: pEqn.H:38
U
Definition: pEqn.H:83
volVectorField F(fluid.F())
surfaceScalarField alphaf1("alphaf1", fvc::interpolate(alpha1))
surfaceScalarField & phi2
multiphaseSystem & fluid
Definition: createFields.H:10
p
Definition: pEqn.H:50
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
phiHbyA
Definition: pEqn.H:20
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
surfaceScalarField & alphaPhi2
const dictionary & pimple
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
surfaceScalarField & phi1
tmp< surfaceScalarField > interpolate(const RhoType &rho)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
forAll(phases, phasei)
Definition: pEqn.H:5
fv::options & fvOptions
surfaceScalarField phir(phic *interface.nHatf())
rho1
Definition: pEqn.H:114
const surfaceScalarField & ghf
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
surfaceScalarField snGradAlpha1(fvc::snGrad(alpha1)*mesh.magSf())
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar pos(const dimensionedScalar &ds)
tmp< surfaceScalarField > DbyA1(fvc::interpolate(rAU1 *(D+phase1.turbulence().pPrime())))
dynamicFvMesh & mesh
volScalarField Kd(fluid.Kd())
rhoEqn solve()
IOMRFZoneList & MRF
rho2
Definition: pEqn.H:115
surfaceScalarField ghSnGradRho(ghf *fvc::snGrad(rho)*mesh.magSf())
surfaceScalarField alpharAUf1(fvc::interpolate(max(alpha1, phase1.residualAlpha())*rAU1))
volVectorField & U1
surfaceScalarField alpharAUf2(fvc::interpolate(max(alpha2, phase2.residualAlpha())*rAU2))
alpha2
Definition: alphaEqn.H:112
phaseModel & phase1
p_rgh
Definition: pEqn.H:120
const dimensionedVector & g
setSnGrad< fixedFluxPressureFvPatchScalarField >(p_rgh.boundaryFieldRef(),(phiHbyA.boundaryField()-MRF.relative(phib))/(mesh.magSf().boundaryField()*rAUf.boundaryField()))
volScalarField & alpha1
dimensionedScalar pMin("pMin", dimPressure, fluid)
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
volScalarField rAU2(IOobject::groupName("rAU", phase2.name()), 1.0/(U2Eqn.A()+max(phase2.residualAlpha()-alpha2, scalar(0))*rho2/runTime.deltaT()))
label patchi
tmp< surfaceScalarField > phiF1
Definition: pEqn.H:37
psi1
Definition: TEqns.H:34
volVectorField & U2
const volScalarField & gh
volScalarField p_rgh_0(p_rgh)
tmp< surfaceScalarField > DbyA2(fvc::interpolate(rAU2 *(D+phase2.turbulence().pPrime())))
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
dimensioned< scalar > mag(const dimensioned< Type > &)
phi
Definition: pEqn.H:18
surfaceScalarField alphaf2("alphaf2", scalar(1)-alphaf1)
phaseModel & phase2
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
rho
Definition: pEqn.H:1
void deleteDemandDrivenData(DataPtr &dataPtr)
surfaceScalarField & alphaPhi1
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45