pEqn.H
Go to the documentation of this file.
2 surfaceScalarField alphaf2("alphaf2", scalar(1) - alphaf1);
3 
5 (
6  "alphaRhof10",
8  (
9  max(alpha1.oldTime(), phase1.residualAlpha())
10  *rho1.oldTime()
11  )
12 );
13 
15 (
16  "alphaRhof20",
18  (
19  max(alpha2.oldTime(), phase2.residualAlpha())
20  *rho2.oldTime()
21  )
22 );
23 
24 // Drag coefficient
25 surfaceScalarField Kdf("Kdf", fluid.Kdf());
26 
27 // Virtual-mass coefficient
28 surfaceScalarField Vmf("Vmf", fluid.Vmf());
29 
31 (
32  IOobject::groupName("rAUf", phase1.name()),
33  1.0
34  /(
35  (alphaRhof10 + Vmf)/runTime.deltaT()
36  + fvc::interpolate(U1Eqn.A())
37  + Kdf
38  )
39 );
40 
42 (
43  IOobject::groupName("rAUf", phase2.name()),
44  1.0
45  /(
46  (alphaRhof20 + Vmf)/runTime.deltaT()
47  + fvc::interpolate(U2Eqn.A())
48  + Kdf
49  )
50 );
51 
52 
53 // Turbulent dispersion, particle-pressure, lift and wall-lubrication forces
54 tmp<surfaceScalarField> Ff1;
55 tmp<surfaceScalarField> Ff2;
56 {
57  // Turbulent-dispersion diffusivity
58  volScalarField D(fluid.D());
59 
60  // Phase-1 turbulent dispersion and particle-pressure diffusivity
62  (
63  fvc::interpolate(D + phase1.turbulence().pPrime())
64  );
65 
66  // Phase-2 turbulent dispersion and particle-pressure diffusivity
68  (
69  fvc::interpolate(D + phase2.turbulence().pPrime())
70  );
71 
72  // Cache the net diffusivity for implicit diffusion treatment in the
73  // phase-fraction equation
74  if (implicitPhasePressure)
75  {
76  fluid.pPrimeByA() = rAUf1*Df1 + rAUf2*Df2;
77  }
78 
79  // Lift and wall-lubrication forces
81 
82  // Phase-fraction face-gradient
84 
85  // Phase-1 dispersion, lift and wall-lubrication force
86  Ff1 = Df1*snGradAlpha1 + Ff;
87 
88  // Phase-2 dispersion, lift and wall-lubrication force
89  Ff2 = -Df2*snGradAlpha1 - Ff;
90 }
91 
92 
93 while (pimple.correct())
94 {
95  // Update continuity errors due to temperature changes
96  #include "correctContErrs.H"
97 
100 
101  // Correct fixed-flux BCs to be consistent with the velocity BCs
102  MRF.correctBoundaryFlux(U1, phi1);
103  MRF.correctBoundaryFlux(U2, phi2);
104 
106  (
107  IOobject::groupName("alpharAUf", phase1.name()),
108  max(alphaf1, phase1.residualAlpha())*rAUf1
109  );
110 
112  (
113  IOobject::groupName("alpharAUf", phase2.name()),
114  max(alphaf2, phase2.residualAlpha())*rAUf2
115  );
116 
117  volScalarField rho("rho", fluid.rho());
119  (
120  "ghSnGradRho",
121  ghf*fvc::snGrad(rho)*mesh.magSf()
122  );
123 
124  // Phase-1 buoyancy flux
125  surfaceScalarField phig1
126  (
127  IOobject::groupName("phig", phase1.name()),
128  alpharAUf1
129  *(
131  - alphaf2*(rhof1 - rhof2)*(g & mesh.Sf())
132  )
133  );
134 
135  // Phase-2 buoyancy flux
136  surfaceScalarField phig2
137  (
138  IOobject::groupName("phig", phase2.name()),
139  alpharAUf2
140  *(
142  - alphaf1*(rhof2 - rhof1)*(g & mesh.Sf())
143  )
144  );
145 
146 
147  // Phase-1 predicted flux
148  surfaceScalarField phiHbyA1
149  (
150  IOobject::groupName("phiHbyA", phase1.name()),
151  phi1
152  );
153 
154  phiHbyA1 =
155  rAUf1
156  *(
157  (alphaRhof10 + Vmf)
158  *MRF.absolute(phi1.oldTime())/runTime.deltaT()
159  + (fvc::interpolate(U1Eqn.H()) & mesh.Sf())
160  + Vmf*ddtPhi2
161  + Kdf*MRF.absolute(phi2)
162  - Ff1()
163  );
164 
165  // Phase-2 predicted flux
166  surfaceScalarField phiHbyA2
167  (
168  IOobject::groupName("phiHbyA", phase2.name()),
169  phi2
170  );
171 
172  phiHbyA2 =
173  rAUf2
174  *(
175  (alphaRhof20 + Vmf)
176  *MRF.absolute(phi2.oldTime())/runTime.deltaT()
177  + (fvc::interpolate(U2Eqn.H()) & mesh.Sf())
178  + Vmf*ddtPhi1
179  + Kdf*MRF.absolute(phi1)
180  - Ff2()
181  );
182 
183 
185  (
186  "phiHbyA",
187  alphaf1*(phiHbyA1 - phig1) + alphaf2*(phiHbyA2 - phig2)
188  );
189  MRF.makeRelative(phiHbyA);
190 
191  phiHbyA1 -= phig1;
192  phiHbyA2 -= phig2;
193 
195  (
196  "rAUf",
198  );
199 
200  // Update the fixedFluxPressure BCs to ensure flux consistency
202  (
203  p_rgh.boundaryField(),
204  (
205  phiHbyA.boundaryField()
206  - (
207  alphaf1.boundaryField()*phi1.boundaryField()
208  + alphaf2.boundaryField()*phi2.boundaryField()
209  )
210  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
211  );
212 
213  tmp<fvScalarMatrix> pEqnComp1;
214  tmp<fvScalarMatrix> pEqnComp2;
215 
216  if (pimple.transonic())
217  {
218  surfaceScalarField phid1
219  (
220  IOobject::groupName("phid", phase1.name()),
222  );
223  surfaceScalarField phid2
224  (
225  IOobject::groupName("phid", phase2.name()),
227  );
228 
229  pEqnComp1 =
230  (
231  contErr1
233  )/rho1
235  (
237  + fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
238  );
239  deleteDemandDrivenData(pEqnComp1().faceFluxCorrectionPtr());
240  pEqnComp1().relax();
241 
242  pEqnComp2 =
243  (
244  contErr2
246  )/rho2
248  (
250  + fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
251  );
252  deleteDemandDrivenData(pEqnComp2().faceFluxCorrectionPtr());
253  pEqnComp2().relax();
254  }
255  else
256  {
257  pEqnComp1 =
258  (
259  contErr1
261  )/rho1
263 
264  pEqnComp2 =
265  (
266  contErr2
268  )/rho2
270  }
271 
272  // Cache p prior to solve for density update
273  volScalarField p_rgh_0("p_rgh_0", p_rgh);
274 
275  while (pimple.correctNonOrthogonal())
276  {
277  fvScalarMatrix pEqnIncomp
278  (
281  );
282 
283  solve
284  (
285  pEqnComp1() + pEqnComp2() + pEqnIncomp,
286  mesh.solver(p_rgh.select(pimple.finalInnerIter()))
287  );
288 
289  if (pimple.finalNonOrthogonalIter())
290  {
291  surfaceScalarField mSfGradp("mSfGradp", pEqnIncomp.flux()/rAUf);
292 
293  phi = phiHbyA + pEqnIncomp.flux();
294 
295  surfaceScalarField phi1s
296  (
297  phiHbyA1
298  + alpharAUf1*mSfGradp
299  - rAUf1*Kdf*MRF.absolute(phi2)
300  );
301 
302  surfaceScalarField phi2s
303  (
304  phiHbyA2
305  + alpharAUf2*mSfGradp
306  - rAUf2*Kdf*MRF.absolute(phi1)
307  );
308 
310  (
311  ((phi2s + rAUf2*Kdf*phi1s) - (phi1s + rAUf1*Kdf*phi2s))
312  /(1.0 - rAUf1*rAUf2*sqr(Kdf))
313  );
314 
315  phi1 = phi - alphaf2*phir;
316  phi2 = phi + alphaf1*phir;
317 
318  U1 = fvc::reconstruct(MRF.absolute(phi1));
319  U1.correctBoundaryConditions();
320  fvOptions.correct(U1);
321 
322  U2 = fvc::reconstruct(MRF.absolute(phi2));
323  U2.correctBoundaryConditions();
324  fvOptions.correct(U2);
325 
326  U = fluid.U();
327 
328  fluid.dgdt() =
329  (
330  alpha1*(pEqnComp2 & p_rgh)
331  - alpha2*(pEqnComp1 & p_rgh)
332  );
333  }
334  }
335 
336  // Update and limit the static pressure
337  p = max(p_rgh + rho*gh, pMin);
338 
339  // Limit p_rgh
340  p_rgh = p - rho*gh;
341 
342  // Update densities from change in p_rgh
343  rho1 += psi1*(p_rgh - p_rgh_0);
344  rho2 += psi2*(p_rgh - p_rgh_0);
345 
346  // Correct p_rgh for consistency with p and the updated densities
347  rho = fluid.rho();
348  p_rgh = p - rho*gh;
349  p_rgh.correctBoundaryConditions();
350 }
351 
352 // Update the phase kinetic energies
353 K1 = 0.5*magSqr(U1);
354 K2 = 0.5*magSqr(U2);
355 
356 // Update the pressure time-derivative if required
357 if (thermo1.dpdt() || thermo2.dpdt())
358 {
359  dpdt = fvc::ddt(p);
360 }
K1
K1
Definition: pEqn.H:403
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
contErr2
contErr2
Definition: correctContErrs.H:5
rAUf2
surfaceScalarField rAUf2(IOobject::groupName("rAUf", phase2.name()), 1.0/((alphaRhof20+Vmf)/runTime.deltaT()+fvc::interpolate(U2Eqn.A())+Kdf))
thermo2
rhoThermo & thermo2
Definition: createFields.H:40
p
p
Definition: pEqn.H:62
Ff
surfaceScalarField Ff(fluid.Ff())
ghf
const surfaceScalarField & ghf
Definition: setRegionFluidFields.H:35
phi1
surfaceScalarField & phi1
Definition: createFields.H:19
snGradAlpha1
surfaceScalarField snGradAlpha1(fvc::snGrad(alpha1) *mesh.magSf())
contErr1
contErr1
Definition: correctContErrs.H:1
fluid
multiphaseSystem & fluid
Definition: createFields.H:10
correctContErrs.H
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
phir
surfaceScalarField phir(phic *interface.nHatf())
rAUf1
surfaceScalarField rAUf1(IOobject::groupName("rAUf", phase1.name()), 1.0/((alphaRhof10+Vmf)/runTime.deltaT()+fvc::interpolate(U1Eqn.A())+Kdf))
dpdt
volScalarField & dpdt
Definition: setRegionFluidFields.H:12
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
phi2
surfaceScalarField & phi2
Definition: createFields.H:24
phi
phi
Definition: pEqn.H:20
U
U
Definition: pEqn.H:46
K2
K2
Definition: pEqn.H:404
Df1
surfaceScalarField Df1(fvc::interpolate(D+phase1.turbulence().pPrime()))
Vmf
surfaceScalarField Vmf("Vmf", fluid.Vmf())
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
ddtPhi2
ddtPhi2
Definition: DDtU.H:6
Ff2
tmp< surfaceScalarField > Ff2
Definition: pEqn.H:55
Df2
surfaceScalarField Df2(fvc::interpolate(D+phase2.turbulence().pPrime()))
psi2
psi2
Definition: TEqns.H:35
thermo1
rhoThermo & thermo1
Definition: createFields.H:39
ddtPhi1
ddtPhi1
Definition: DDtU.H:1
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Ff1
tmp< surfaceScalarField > Ff1
Definition: pEqn.H:54
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
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::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
ghSnGradRho
surfaceScalarField ghSnGradRho(ghf *fvc::snGrad(rho) *mesh.magSf())
phase2
phaseModel & phase2
Definition: createFields.H:13
alphaf2
surfaceScalarField alphaf2("alphaf2", scalar(1) - alphaf1)
alphaRhof20
surfaceScalarField alphaRhof20("alphaRhof20", fvc::interpolate(max(alpha2.oldTime(), phase2.residualAlpha()) *rho2.oldTime()))
alphaPhi2
surfaceScalarField & alphaPhi2
Definition: createFields.H:25
alphaRhof10
surfaceScalarField alphaRhof10("alphaRhof10", fvc::interpolate(max(alpha1.oldTime(), phase1.residualAlpha()) *rho1.oldTime()))
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
Kdf
surfaceScalarField Kdf("Kdf", fluid.Kdf())
p_rgh
volScalarField & p_rgh
Definition: setRegionFluidFields.H:31
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:16
phase1
phaseModel & phase1
Definition: createFields.H:12
rAUf
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))