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 phase diffusivities for implicit treatment in the
73  // phase-fraction equation
74  if (implicitPhasePressure)
75  {
76  phase1.DbyA(rAUf1*Df1);
77  phase2.DbyA(rAUf2*Df2);
78  }
79 
80  // Lift and wall-lubrication forces
82 
83  // Phase-fraction face-gradient
85 
86  // Phase-1 dispersion, lift and wall-lubrication force
87  Ff1 = Df1*snGradAlpha1 + Ff;
88 
89  // Phase-2 dispersion, lift and wall-lubrication force
90  Ff2 = -Df2*snGradAlpha1 - Ff;
91 }
92 
93 
94 while (pimple.correct())
95 {
96  // Update continuity errors due to temperature changes
97  fluid.correct();
98 
101 
102  // Correct fixed-flux BCs to be consistent with the velocity BCs
103  MRF.correctBoundaryFlux(U1, phi1);
104  MRF.correctBoundaryFlux(U2, phi2);
105 
107  (
108  IOobject::groupName("alpharAUf", phase1.name()),
109  max(alphaf1, phase1.residualAlpha())*rAUf1
110  );
111 
113  (
114  IOobject::groupName("alpharAUf", phase2.name()),
115  max(alphaf2, phase2.residualAlpha())*rAUf2
116  );
117 
118  volScalarField rho("rho", fluid.rho());
120  (
121  "ghSnGradRho",
122  ghf*fvc::snGrad(rho)*mesh.magSf()
123  );
124 
125  // Phase-1 buoyancy flux
126  surfaceScalarField phig1
127  (
128  IOobject::groupName("phig", phase1.name()),
129  alpharAUf1
130  *(
132  - alphaf2*(rhof1 - rhof2)*(g & mesh.Sf())
133  )
134  );
135 
136  // Phase-2 buoyancy flux
137  surfaceScalarField phig2
138  (
139  IOobject::groupName("phig", phase2.name()),
140  alpharAUf2
141  *(
143  - alphaf1*(rhof2 - rhof1)*(g & mesh.Sf())
144  )
145  );
146 
147 
148  // Phase-1 predicted flux
149  surfaceScalarField phiHbyA1
150  (
151  IOobject::groupName("phiHbyA", phase1.name()),
152  phi1
153  );
154 
155  phiHbyA1 =
156  rAUf1
157  *(
158  (alphaRhof10 + Vmf)
159  *MRF.absolute(phi1.oldTime())/runTime.deltaT()
160  + (fvc::interpolate(U1Eqn.H()) & mesh.Sf())
161  + Vmf*ddtPhi2
162  + Kdf*MRF.absolute(phi2)
163  - Ff1()
164  );
165 
166  // Phase-2 predicted flux
167  surfaceScalarField phiHbyA2
168  (
169  IOobject::groupName("phiHbyA", phase2.name()),
170  phi2
171  );
172 
173  phiHbyA2 =
174  rAUf2
175  *(
176  (alphaRhof20 + Vmf)
177  *MRF.absolute(phi2.oldTime())/runTime.deltaT()
178  + (fvc::interpolate(U2Eqn.H()) & mesh.Sf())
179  + Vmf*ddtPhi1
180  + Kdf*MRF.absolute(phi1)
181  - Ff2()
182  );
183 
184 
186  (
187  "phiHbyA",
188  alphaf1*(phiHbyA1 - phig1) + alphaf2*(phiHbyA2 - phig2)
189  );
190  MRF.makeRelative(phiHbyA);
191 
192  phiHbyA1 -= phig1;
193  phiHbyA2 -= phig2;
194 
196  (
197  "rAUf",
199  );
200 
201  // Update the fixedFluxPressure BCs to ensure flux consistency
203  (
204  p_rgh.boundaryField(),
205  (
206  phiHbyA.boundaryField()
207  - (
208  alphaf1.boundaryField()*phi1.boundaryField()
209  + alphaf2.boundaryField()*phi2.boundaryField()
210  )
211  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
212  );
213 
214  tmp<fvScalarMatrix> pEqnComp1;
215  tmp<fvScalarMatrix> pEqnComp2;
216 
217  if (pimple.transonic())
218  {
219  surfaceScalarField phid1
220  (
221  IOobject::groupName("phid", phase1.name()),
223  );
224  surfaceScalarField phid2
225  (
226  IOobject::groupName("phid", phase2.name()),
228  );
229 
230  if (phase1.compressible())
231  {
232  pEqnComp1 =
233  (
234  phase1.continuityError()
236  )/rho1
238  (
240  + fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
241  );
242  deleteDemandDrivenData(pEqnComp1().faceFluxCorrectionPtr());
243  pEqnComp1().relax();
244  }
245 
246  if (phase2.compressible())
247  {
248  pEqnComp2 =
249  (
250  phase2.continuityError()
252  )/rho2
254  (
256  + fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
257  );
258  deleteDemandDrivenData(pEqnComp2().faceFluxCorrectionPtr());
259  pEqnComp2().relax();
260  }
261  }
262  else
263  {
264  if (phase1.compressible())
265  {
266  pEqnComp1 =
267  (
268  phase1.continuityError()
270  )/rho1
272  }
273 
274  if (phase2.compressible())
275  {
276  pEqnComp2 =
277  (
278  phase2.continuityError()
280  )/rho2
282  }
283  }
284 
285  if (fluid.transfersMass())
286  {
287  if (pEqnComp1.valid())
288  {
289  pEqnComp1() -= fluid.dmdt()/rho1;
290  }
291  else
292  {
293  pEqnComp1 = fvm::Su(-fluid.dmdt()/rho1, p_rgh);
294  }
295 
296  if (pEqnComp2.valid())
297  {
298  pEqnComp2() += fluid.dmdt()/rho2;
299  }
300  else
301  {
302  pEqnComp2 = fvm::Su(fluid.dmdt()/rho2, p_rgh);
303  }
304  }
305 
306  // Cache p prior to solve for density update
307  volScalarField p_rgh_0("p_rgh_0", p_rgh);
308 
309  while (pimple.correctNonOrthogonal())
310  {
311  fvScalarMatrix pEqnIncomp
312  (
315  );
316 
317  {
318  fvScalarMatrix pEqn(pEqnIncomp);
319 
320  if (pEqnComp1.valid())
321  {
322  pEqn += pEqnComp1();
323  }
324 
325  if (pEqnComp2.valid())
326  {
327  pEqn += pEqnComp2();
328  }
329 
330  solve
331  (
332  pEqn,
333  mesh.solver(p_rgh.select(pimple.finalInnerIter()))
334  );
335  }
336 
337  if (pimple.finalNonOrthogonalIter())
338  {
339  surfaceScalarField mSfGradp("mSfGradp", pEqnIncomp.flux()/rAUf);
340 
341  phi = phiHbyA + pEqnIncomp.flux();
342 
343  surfaceScalarField phi1s
344  (
345  phiHbyA1
346  + alpharAUf1*mSfGradp
347  - rAUf1*Kdf*MRF.absolute(phi2)
348  );
349 
350  surfaceScalarField phi2s
351  (
352  phiHbyA2
353  + alpharAUf2*mSfGradp
354  - rAUf2*Kdf*MRF.absolute(phi1)
355  );
356 
358  (
359  ((phi2s + rAUf2*Kdf*phi1s) - (phi1s + rAUf1*Kdf*phi2s))
360  /(1.0 - rAUf1*rAUf2*sqr(Kdf))
361  );
362 
363  phi1 = phi - alphaf2*phir;
364  phi2 = phi + alphaf1*phir;
365 
366  U1 = fvc::reconstruct(MRF.absolute(phi1));
367  U1.correctBoundaryConditions();
368  fvOptions.correct(U1);
369 
370  U2 = fvc::reconstruct(MRF.absolute(phi2));
371  U2.correctBoundaryConditions();
372  fvOptions.correct(U2);
373 
374  // Set the phase dilatation rates
375  if (pEqnComp1.valid())
376  {
377  phase1.divU(-pEqnComp1 & p_rgh);
378  }
379  if (pEqnComp2.valid())
380  {
381  phase2.divU(-pEqnComp2 & p_rgh);
382  }
383  }
384  }
385 
386  Info<< "min(p) = " << min(p_rgh + rho*gh).value() << endl;
387  if (min(p_rgh + rho*gh) < pMin)
388  {
389  Info<< "Clipping p" << endl;
390  }
391 
392  // Update and limit the static pressure
393  p = max(p_rgh + rho*gh, pMin);
394 
395  // Limit p_rgh
396  p_rgh = p - rho*gh;
397 
398  // Update densities from change in p_rgh
399  rho1 += psi1*(p_rgh - p_rgh_0);
400  rho2 += psi2*(p_rgh - p_rgh_0);
401 
402  // Correct p_rgh for consistency with p and the updated densities
403  rho = fluid.rho();
404  p_rgh = p - rho*gh;
405  p_rgh.correctBoundaryConditions();
406 }
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
rAUf2
surfaceScalarField rAUf2(IOobject::groupName("rAUf", phase2.name()), 1.0/((alphaRhof20+Vmf)/runTime.deltaT()+fvc::interpolate(U2Eqn.A())+Kdf))
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())
fluid
multiphaseSystem & fluid
Definition: createFields.H:10
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::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::correction
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
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))
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
phi2
surfaceScalarField & phi2
Definition: createFields.H:24
phi
phi
Definition: pEqn.H:20
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
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
ddtPhi2
ddtPhi2
Definition: DDtU.H:6
Foam::Info
messageStream Info
Ff2
tmp< surfaceScalarField > Ff2
Definition: pEqn.H:55
Df2
surfaceScalarField Df2(fvc::interpolate(D+phase2.turbulence().pPrime()))
psi2
psi2
Definition: TEqns.H:35
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())
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
p_rgh
volScalarField & p_rgh
Definition: setRegionFluidFields.H:31
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:16
phase1
phaseModel & phase1
Definition: createFields.H:12
rAUf
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))