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 // Turbulent diffusion, particle-pressure, lift and wall-lubrication fluxes
35 tmp<surfaceScalarField> phiF1;
36 tmp<surfaceScalarField> phiF2;
37 {
38  // Turbulent-dispersion diffusivity
39  volScalarField D(fluid.D());
40 
41  // Phase-1 turbulent dispersion and particle-pressure flux
43  (
45  (
46  rAU1*(D + phase1.turbulence().pPrime())
47  )
48  );
49 
50  // Phase-2 turbulent dispersion and particle-pressure flux
52  (
54  (
55  rAU2*(D + phase2.turbulence().pPrime())
56  )
57  );
58 
59  // Cache the net diffusivity for implicit diffusion treatment in the
60  // phase-fraction equation
61  if (implicitPhasePressure)
62  {
63  fluid.pPrimeByA() = Df1 + Df2;
64  }
65 
66  // Lift and wall-lubrication forces
67  volVectorField F(fluid.F());
68 
69  // Phase-fraction face-gradient
71 
72  // Phase-1 dispersion, lift and wall-lubrication flux
73  phiF1 = Df1*snGradAlpha1 + fvc::flux(rAU1*F);
74 
75  // Phase-1 dispersion, lift and wall-lubrication flux
76  phiF2 = - Df2*snGradAlpha1 - fvc::flux(rAU2*F);
77 }
78 
79 
80 // --- Pressure corrector loop
81 while (pimple.correct())
82 {
83  // Update continuity errors due to temperature changes
84  #include "correctContErrs.H"
85 
86  volScalarField rho("rho", fluid.rho());
87 
88  // Correct p_rgh for consistency with p and the updated densities
89  p_rgh = p - rho*gh;
90 
91  // Correct fixed-flux BCs to be consistent with the velocity BCs
92  MRF.correctBoundaryFlux(U1, phi1);
93  MRF.correctBoundaryFlux(U2, phi2);
94 
95  volVectorField HbyA1
96  (
97  IOobject::groupName("HbyA", phase1.name()),
98  U1
99  );
100  HbyA1 =
101  rAU1
102  *(
103  U1Eqn.H()
104  + max(phase1.residualAlpha() - alpha1, scalar(0))
105  *rho1*U1.oldTime()/runTime.deltaT()
106  );
107 
108  volVectorField HbyA2
109  (
110  IOobject::groupName("HbyA", phase2.name()),
111  U2
112  );
113  HbyA2 =
114  rAU2
115  *(
116  U2Eqn.H()
117  + max(phase2.residualAlpha() - alpha2, scalar(0))
118  *rho2*U2.oldTime()/runTime.deltaT()
119  );
120 
122  (
123  "ghSnGradRho",
124  ghf*fvc::snGrad(rho)*mesh.magSf()
125  );
126 
127  surfaceScalarField phig1
128  (
129  alpharAUf1
130  *(
132  - alphaf2*fvc::interpolate(rho1 - rho2)*(g & mesh.Sf())
133  )
134  );
135 
136  surfaceScalarField phig2
137  (
138  alpharAUf2
139  *(
141  - alphaf1*fvc::interpolate(rho2 - rho1)*(g & mesh.Sf())
142  )
143  );
144 
145 
146  // ddtPhiCorr filter -- only apply in pure(ish) phases
148  surfaceScalarField phiCorrCoeff1(pos(alphaf1Bar - 0.99));
149  surfaceScalarField phiCorrCoeff2(pos(0.01 - alphaf1Bar));
150 
151  {
152  surfaceScalarField::Boundary& phiCorrCoeff1Bf =
153  phiCorrCoeff1.boundaryFieldRef();
154 
155  surfaceScalarField::Boundary& phiCorrCoeff2Bf =
156  phiCorrCoeff2.boundaryFieldRef();
157 
158  forAll(mesh.boundary(), patchi)
159  {
160  // Set ddtPhiCorr to 0 on non-coupled boundaries
161  if
162  (
163  !mesh.boundary()[patchi].coupled()
164  || isA<cyclicAMIFvPatch>(mesh.boundary()[patchi])
165  )
166  {
167  phiCorrCoeff1Bf[patchi] = 0;
168  phiCorrCoeff2Bf[patchi] = 0;
169  }
170  }
171  }
172 
173  // Phase-1 predicted flux
174  surfaceScalarField phiHbyA1
175  (
176  IOobject::groupName("phiHbyA", phase1.name()),
177  fvc::flux(HbyA1)
178  + phiCorrCoeff1*fvc::interpolate(alpha1.oldTime()*rho1.oldTime()*rAU1)
179  *(
180  MRF.absolute(phi1.oldTime())
181  - fvc::flux(U1.oldTime())
182  )/runTime.deltaT()
183  - phiF1()
184  - phig1
185  );
186 
187  // Phase-2 predicted flux
188  surfaceScalarField phiHbyA2
189  (
190  IOobject::groupName("phiHbyA", phase2.name()),
191  fvc::flux(HbyA2)
192  + phiCorrCoeff2*fvc::interpolate(alpha2.oldTime()*rho2.oldTime()*rAU2)
193  *(
194  MRF.absolute(phi2.oldTime())
195  - fvc::flux(U2.oldTime())
196  )/runTime.deltaT()
197  - phiF2()
198  - phig2
199  );
200 
201  // Face-drag coefficients
204 
205  // Construct the mean predicted flux
206  // including explicit drag contributions based on absolute fluxes
208  (
209  "phiHbyA",
210  alphaf1*(phiHbyA1 + rAUKd1*MRF.absolute(phi2))
211  + alphaf2*(phiHbyA2 + rAUKd2*MRF.absolute(phi1))
212  );
213  MRF.makeRelative(phiHbyA);
214 
215  // Construct pressure "diffusivity"
217  (
218  "rAUf",
220  );
221 
222  // Update the fixedFluxPressure BCs to ensure flux consistency
224  (
225  p_rgh.boundaryFieldRef(),
226  (
227  phiHbyA.boundaryField()
228  - (
229  alphaf1.boundaryField()*phi1.boundaryField()
230  + alphaf2.boundaryField()*phi2.boundaryField()
231  )
232  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
233  );
234 
235  tmp<fvScalarMatrix> pEqnComp1;
236  tmp<fvScalarMatrix> pEqnComp2;
237 
238  // Construct the compressibility parts of the pressure equation
239  if (pimple.transonic())
240  {
241  surfaceScalarField phid1
242  (
243  IOobject::groupName("phid", phase1.name()),
245  );
246  surfaceScalarField phid2
247  (
248  IOobject::groupName("phid", phase2.name()),
250  );
251 
252  pEqnComp1 =
253  (
254  contErr1
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  deleteDemandDrivenData(pEqnComp1.ref().faceFluxCorrectionPtr());
266  pEqnComp1.ref().relax();
267 
268  pEqnComp2 =
269  (
270  contErr2
272  )/rho2
273  + correction
274  (
275  (alpha2/rho2)*
276  (
277  psi2*fvm::ddt(p_rgh)
278  + fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
279  )
280  );
281  deleteDemandDrivenData(pEqnComp2.ref().faceFluxCorrectionPtr());
282  pEqnComp2.ref().relax();
283  }
284  else
285  {
286  pEqnComp1 =
287  (
288  contErr1
290  )/rho1
291  + (alpha1*psi1/rho1)*correction(fvm::ddt(p_rgh));
292 
293  pEqnComp2 =
294  (
295  contErr2
297  )/rho2
298  + (alpha2*psi2/rho2)*correction(fvm::ddt(p_rgh));
299  }
300 
301  // Cache p prior to solve for density update
302  volScalarField p_rgh_0(p_rgh);
303 
304  // Iterate over the pressure equation to correct for non-orthogonality
305  while (pimple.correctNonOrthogonal())
306  {
307  // Construct the transport part of the pressure equation
308  fvScalarMatrix pEqnIncomp
309  (
310  fvc::div(phiHbyA)
311  - fvm::laplacian(rAUf, p_rgh)
312  );
313 
314  solve
315  (
316  pEqnComp1() + pEqnComp2() + pEqnIncomp,
317  mesh.solver(p_rgh.select(pimple.finalInnerIter()))
318  );
319 
320  // Correct fluxes and velocities on last non-orthogonal iteration
321  if (pimple.finalNonOrthogonalIter())
322  {
323  phi = phiHbyA + pEqnIncomp.flux();
324 
325  surfaceScalarField mSfGradp("mSfGradp", pEqnIncomp.flux()/rAUf);
326 
327  // Partial-elimination phase-flux corrector
328  {
329  surfaceScalarField phi1s
330  (
331  phiHbyA1 + alpharAUf1*mSfGradp
332  );
333 
334  surfaceScalarField phi2s
335  (
336  phiHbyA2 + alpharAUf2*mSfGradp
337  );
338 
340  (
341  ((phi1s + rAUKd1*phi2s) - (phi2s + rAUKd2*phi1s))
342  /(1 - rAUKd1*rAUKd2)
343  );
344 
345  phi1 = phi + alphaf2*phir;
346  phi2 = phi - alphaf1*phir;
347  }
348 
349  // Compressibility correction for phase-fraction equations
350  fluid.dgdt() =
351  (
352  alpha1*(pEqnComp2 & p_rgh)
353  - alpha2*(pEqnComp1 & p_rgh)
354  );
355 
356  // Optionally relax pressure for velocity correction
357  p_rgh.relax();
358 
359  mSfGradp = pEqnIncomp.flux()/rAUf;
360 
361  // Partial-elimination phase-velocity corrector
362  {
363  volVectorField Us1
364  (
365  HbyA1
366  + fvc::reconstruct(alpharAUf1*mSfGradp - phiF1() - phig1)
367  );
368 
369  volVectorField Us2
370  (
371  HbyA2
372  + fvc::reconstruct(alpharAUf2*mSfGradp - phiF2() - phig2)
373  );
374 
375  volScalarField D1(rAU1*Kd);
376  volScalarField D2(rAU2*Kd);
377 
378  U = alpha1*(Us1 + D1*U2) + alpha2*(Us2 + D2*U1);
379  volVectorField Ur(((1 - D2)*Us1 - (1 - D1)*Us2)/(1 - D1*D2));
380 
381  U1 = U + alpha2*Ur;
382  U1.correctBoundaryConditions();
383  fvOptions.correct(U1);
384 
385  U2 = U - alpha1*Ur;
386  U2.correctBoundaryConditions();
387  fvOptions.correct(U2);
388 
389  U = fluid.U();
390  }
391  }
392  }
393 
394  // Update and limit the static pressure
395  p = max(p_rgh + rho*gh, pMin);
396 
397  // Limit p_rgh
398  p_rgh = p - rho*gh;
399 
400  // Update densities from change in p_rgh
401  rho1 += psi1*(p_rgh - p_rgh_0);
402  rho2 += psi2*(p_rgh - p_rgh_0);
403 
404  // Correct p_rgh for consistency with p and the updated densities
405  rho = fluid.rho();
406  p_rgh = p - rho*gh;
407  p_rgh.correctBoundaryConditions();
408 }
409 
410 // Update the phase kinetic energies
411 K1 = 0.5*magSqr(U1);
412 K2 = 0.5*magSqr(U2);
413 
414 // Update the pressure time-derivative if required
415 if (thermo1.dpdt() || thermo2.dpdt())
416 {
417  dpdt = fvc::ddt(p);
418 }
rhoThermo & thermo2
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.
contErr1
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 > &)
surfaceScalarField Df2(fvc::interpolate(D+phase2.turbulence().pPrime()))
contErr2
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 Df1(fvc::interpolate(D+phase1.turbulence().pPrime()))
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
rhoThermo & thermo1
dimensionedScalar pos(const dimensionedScalar &ds)
dynamicFvMesh & mesh
volScalarField Kd(fluid.Kd())
rhoEqn solve()
IOMRFZoneList & MRF
rho2
Definition: pEqn.H:115
surfaceScalarField ghSnGradRho(ghf *fvc::snGrad(rho)*mesh.magSf())
volScalarField & dpdt
surfaceScalarField alpharAUf1(fvc::interpolate(max(alpha1, phase1.residualAlpha())*rAU1))
K2
Definition: pEqn.H:412
volVectorField & U1
surfaceScalarField alpharAUf2(fvc::interpolate(max(alpha2, phase2.residualAlpha())*rAU2))
alpha2
Definition: alphaEqn.H:112
phaseModel & phase1
p_rgh
Definition: pEqn.H:120
dimensioned< scalar > magSqr(const dimensioned< Type > &)
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)
K1
Definition: pEqn.H:411
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< 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
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