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