pEqn.H
Go to the documentation of this file.
1// Face volume fractions
2PtrList<surfaceScalarField> alphafs(phases.size());
4{
5 phaseModel& phase = phases[phasei];
6 const volScalarField& alpha = phase;
7
8 alphafs.set(phasei, fvc::interpolate(alpha).ptr());
9 alphafs[phasei].rename("pEqn" + alphafs[phasei].name());
10}
11
12// Diagonal coefficients
13PtrList<volScalarField> rAUs(phases.size());
14forAll(fluid.movingPhases(), movingPhasei)
15{
16 phaseModel& phase = fluid.movingPhases()[movingPhasei];
17 const volScalarField& alpha = phase;
18
19 rAUs.set
20 (
21 phase.index(),
22 new volScalarField
23 (
24 IOobject::groupName("rAU", phase.name()),
25 1.0
26 /(
27 UEqns[phase.index()].A()
28 + byDt(max(phase.residualAlpha() - alpha, scalar(0))*phase.rho())
29 )
30 )
31 );
32}
33fluid.fillFields("rAU", dimTime/dimDensity, rAUs);
34
35// Phase diagonal coefficients
36PtrList<surfaceScalarField> alpharAUfs(phases.size());
38{
39 phaseModel& phase = phases[phasei];
40 const volScalarField& alpha = phase;
41
42 alpharAUfs.set
43 (
44 phasei,
45 (
46 fvc::interpolate(max(alpha, phase.residualAlpha())*rAUs[phasei])
47 ).ptr()
48 );
49}
50
51// Explicit force fluxes
52PtrList<surfaceScalarField> phiFs(fluid.phiFs(rAUs));
53
54// --- Pressure corrector loop
55while (pimple.correct())
56{
57 volScalarField rho("rho", fluid.rho());
58
59 // Correct p_rgh for consistency with p and the updated densities
60 p_rgh = p - rho*gh;
61
62 // Correct fixed-flux BCs to be consistent with the velocity BCs
63 forAll(fluid.movingPhases(), movingPhasei)
64 {
65 phaseModel& phase = fluid.movingPhases()[movingPhasei];
66 MRF.correctBoundaryFlux(phase.U(), phase.phiRef());
67 }
68
69 // Combined buoyancy and force fluxes
70 PtrList<surfaceScalarField> phigFs(phases.size());
71 {
72 const surfaceScalarField ghSnGradRho
73 (
74 "ghSnGradRho",
75 ghf*fvc::snGrad(rho)*mesh.magSf()
76 );
77
79 {
80 phaseModel& phase = phases[phasei];
81
82 phigFs.set
83 (
84 phasei,
85 (
87 *(
89 - (fvc::interpolate(phase.rho() - rho))*(g & mesh.Sf())
90 - fluid.surfaceTension(phase)*mesh.magSf()
91 )
92 ).ptr()
93 );
94
95 if (phiFs.set(phasei))
96 {
97 phigFs[phasei] += phiFs[phasei];
98 }
99 }
100 }
101
102 // Predicted velocities and fluxes for each phase
103 PtrList<volVectorField> HbyAs(phases.size());
104 PtrList<surfaceScalarField> phiHbyAs(phases.size());
105 {
106 // Correction force fluxes
107 PtrList<surfaceScalarField> ddtCorrByAs(fluid.ddtCorrByAs(rAUs));
108
109 forAll(fluid.movingPhases(), movingPhasei)
110 {
111 phaseModel& phase = fluid.movingPhases()[movingPhasei];
112 const volScalarField& alpha = phase;
113
114 HbyAs.set
115 (
116 phase.index(),
117 new volVectorField
118 (
119 IOobject::groupName("HbyA", phase.name()),
120 phase.U()
121 )
122 );
123
124 HbyAs[phase.index()] =
125 rAUs[phase.index()]
126 *(
127 UEqns[phase.index()].H()
128 + byDt
129 (
130 max(phase.residualAlpha() - alpha, scalar(0))
131 *phase.rho()
132 )
133 *phase.U()().oldTime()
134 );
135
136 phiHbyAs.set
137 (
138 phase.index(),
139 new surfaceScalarField
140 (
141 IOobject::groupName("phiHbyA", phase.name()),
142 fvc::flux(HbyAs[phase.index()])
143 - phigFs[phase.index()]
144 - ddtCorrByAs[phase.index()]
145 )
146 );
147 }
148 }
149 fluid.fillFields("HbyA", dimVelocity, HbyAs);
150 fluid.fillFields("phiHbyA", dimForce/dimDensity/dimVelocity, phiHbyAs);
151
152 // Add explicit drag forces and fluxes if not doing partial elimination
153 if (!partialElimination)
154 {
155 PtrList<volVectorField> KdUByAs(fluid.KdUByAs(rAUs));
156 PtrList<surfaceScalarField> phiKdPhis(fluid.phiKdPhis(rAUs));
157
159 {
160 if (KdUByAs.set(phasei))
161 {
162 HbyAs[phasei] -= KdUByAs[phasei];
163 }
164 if (phiKdPhis.set(phasei))
165 {
166 phiHbyAs[phasei] -= phiKdPhis[phasei];
167 }
168 }
169 }
170
171 // Total predicted flux
172 surfaceScalarField phiHbyA
173 (
174 IOobject
175 (
176 "phiHbyA",
177 runTime.timeName(),
178 mesh,
179 IOobject::NO_READ,
180 IOobject::AUTO_WRITE
181 ),
182 mesh,
183 dimensionedScalar("phiHbyA", dimArea*dimVelocity, 0)
184 );
185
187 {
189 }
190
191 // Add explicit drag fluxes if doing partial elimination
192 if (partialElimination)
193 {
194 PtrList<surfaceScalarField> phiKdPhis(fluid.phiKdPhis(rAUs));
195
197 {
198 if (phiKdPhis.set(phasei))
199 {
200 phiHbyA -= alphafs[phasei]*phiKdPhis[phasei];
201 }
202 }
203 }
204
205 MRF.makeRelative(phiHbyA);
206
207 // Construct pressure "diffusivity"
208 surfaceScalarField rAUf
209 (
210 IOobject
211 (
212 "rAUf",
213 runTime.timeName(),
214 mesh
215 ),
216 mesh,
217 dimensionedScalar("rAUf", dimensionSet(-1, 3, 1, 0, 0), 0)
218 );
219
221 {
223 }
224
225 rAUf = mag(rAUf);
226
227 // Update the fixedFluxPressure BCs to ensure flux consistency
228 {
229 surfaceScalarField::Boundary phib(phi.boundaryField());
230 phib = 0;
232 {
233 phaseModel& phase = phases[phasei];
234 phib +=
235 alphafs[phasei].boundaryField()*phase.phi()().boundaryField();
236 }
237
239 (
240 p_rgh.boundaryFieldRef(),
241 (
242 phiHbyA.boundaryField() - phib
243 )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
244 );
245 }
246
247 // Compressible pressure equations
248 PtrList<fvScalarMatrix> pEqnComps(phases.size());
249 PtrList<volScalarField> dmdts(fluid.dmdts());
251 {
252 phaseModel& phase = phases[phasei];
253 const volScalarField& alpha = phase;
254 volScalarField& rho = phase.thermoRef().rho();
255
256 if (phase.compressible())
257 {
258 if (pimple.transonic())
259 {
260 surfaceScalarField phid
261 (
262 IOobject::groupName("phid", phase.name()),
263 fvc::interpolate(phase.thermo().psi())*phase.phi()
264 );
265
266 pEqnComps.set
267 (
268 phasei,
269 (
270 (
271 fvc::ddt(alpha, rho) + fvc::div(phase.alphaRhoPhi())
272 - fvc::Sp
273 (
274 fvc::ddt(alpha) + fvc::div(phase.alphaPhi()),
275 rho
276 )
277 )/rho
278 + correction
279 (
280 (alpha/rho)*
281 (
282 phase.thermo().psi()*fvm::ddt(p_rgh)
283 + fvm::div(phid, p_rgh)
284 - fvm::Sp(fvc::div(phid), p_rgh)
285 )
286 )
287 ).ptr()
288 );
289
290 pEqnComps[phasei].relax();
291 }
292 else
293 {
294 pEqnComps.set
295 (
296 phasei,
297 (
298 (
299 fvc::ddt(alpha, rho) + fvc::div(phase.alphaRhoPhi())
300 - fvc::Sp
301 (
302 (fvc::ddt(alpha) + fvc::div(phase.alphaPhi())),
303 rho
304 )
305 )/rho
306 + (alpha*phase.thermo().psi()/rho)
307 *correction(fvm::ddt(p_rgh))
308 ).ptr()
309 );
310 }
311 }
312
313 // Add option sources
314 if (fvOptions.appliesToField(rho.name()))
315 {
316 tmp<fvScalarMatrix> optEqn = fvOptions(alpha, rho);
317 if (pEqnComps.set(phasei))
318 {
319 pEqnComps[phasei] -= (optEqn&rho)/rho;
320 }
321 else
322 {
323 pEqnComps.set
324 (
325 phasei,
326 fvm::Su(- (optEqn&rho)/rho, p_rgh).ptr()
327 );
328 }
329 }
330
331 // Add mass transfer
332 if (dmdts.set(phasei))
333 {
334 if (pEqnComps.set(phasei))
335 {
336 pEqnComps[phasei] -= dmdts[phasei]/rho;
337 }
338 else
339 {
340 pEqnComps.set
341 (
342 phasei,
343 fvm::Su(- dmdts[phasei]/rho, p_rgh)
344 );
345 }
346 }
347 }
348
349 // Cache p prior to solve for density update
350 volScalarField p_rgh_0(p_rgh);
351
352 // Iterate over the pressure equation to correct for non-orthogonality
353 while (pimple.correctNonOrthogonal())
354 {
355 // Construct the transport part of the pressure equation
356 fvScalarMatrix pEqnIncomp
357 (
358 fvc::div(phiHbyA)
359 - fvm::laplacian(rAUf, p_rgh)
360 );
361
362 {
363 fvScalarMatrix pEqn(pEqnIncomp);
364
366 {
367 if (pEqnComps.set(phasei))
368 {
369 pEqn += pEqnComps[phasei];
370 }
371 }
372
373 solve
374 (
375 pEqn,
376 mesh.solver(p_rgh.select(pimple.finalInnerIter()))
377 );
378 }
379
380 // Correct fluxes and velocities on last non-orthogonal iteration
381 if (pimple.finalNonOrthogonalIter())
382 {
383 phi = phiHbyA + pEqnIncomp.flux();
384
385 surfaceScalarField mSfGradp("mSfGradp", pEqnIncomp.flux()/rAUf);
386
387 forAll(fluid.movingPhases(), movingPhasei)
388 {
389 phaseModel& phase = fluid.movingPhases()[movingPhasei];
390
391 phase.phiRef() =
392 phiHbyAs[phase.index()]
393 + alpharAUfs[phase.index()]*mSfGradp;
394
395 // Set the phase dilatation rates
396 if (pEqnComps.set(phase.index()))
397 {
398 phase.divU(-pEqnComps[phase.index()] & p_rgh);
399 }
400 }
401
402 // Optionally relax pressure for velocity correction
403 p_rgh.relax();
404
405 mSfGradp = pEqnIncomp.flux()/rAUf;
406
407 forAll(fluid.movingPhases(), movingPhasei)
408 {
409 phaseModel& phase = fluid.movingPhases()[movingPhasei];
410
411 phase.URef() =
412 HbyAs[phase.index()]
413 + fvc::reconstruct
414 (
415 alpharAUfs[phase.index()]*mSfGradp
416 - phigFs[phase.index()]
417 );
418 }
419
420 if (partialElimination)
421 {
422 fluid.partialElimination(rAUs);
423 }
424
425 forAll(fluid.movingPhases(), movingPhasei)
426 {
427 phaseModel& phase = fluid.movingPhases()[movingPhasei];
428
429 phase.URef().correctBoundaryConditions();
430 fvOptions.correct(phase.URef());
431 }
432 }
433 }
434
435 // Update and limit the static pressure
436 p = max(p_rgh + rho*gh, pMin);
437
438 // Limit p_rgh
439 p_rgh = p - rho*gh;
440
441 // Update densities from change in p_rgh
443 {
444 phaseModel& phase = phases[phasei];
445 phase.thermoRef().rho() += phase.thermo().psi()*(p_rgh - p_rgh_0);
446 }
447
448 // Correct p_rgh for consistency with p and the updated densities
449 rho = fluid.rho();
450 p_rgh = p - rho*gh;
451 p_rgh.correctBoundaryConditions();
452}
Y[inertIndex] max(0.0)
const uniformDimensionedVectorField & g
volScalarField & p_rgh
fv::options & fvOptions
surfaceScalarField & phi
const surfaceScalarField & ghf
IOMRFZoneList & MRF
const volScalarField & gh
pimpleControl & pimple
twoPhaseSystem & fluid
const dimensionedScalar & pMin
volScalarField & p
phiHbyA
Definition: pcEqn.H:73
surfaceScalarField phid("phid", fvc::interpolate(psi) *(fvc::flux(HbyA)+MRF.zeroFilter(rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho))))
dynamicFvMesh & mesh
engineTime & runTime
Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));if(U.nOldTimes()){ volVectorField *Uold=&U.oldTime();volScalarField *Kold=&K.oldTime();*Kold==0.5 *magSqr(*Uold);while(Uold->nOldTimes()) { Uold=&Uold-> oldTime()
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
PtrList< fvVectorMatrix > UEqns(fluid.phases().size())
volScalarField p_rgh_0(p_rgh)
label phasei
Definition: pEqn.H:27
phib
Definition: pEqn.H:189
surfaceScalarField ghSnGradRho(ghf *fvc::snGrad(rho) *mesh.magSf())
PtrList< surfaceScalarField > phiHbyAs(fluid.phases().size())
PtrList< volVectorField > HbyAs(fluid.phases().size())
setSnGrad< fixedFluxPressureFvPatchScalarField >(p_rgh.boundaryFieldRef(),(phiHbyA.boundaryField() - MRF.relative(phib))/(mesh.magSf().boundaryField() *rAUf.boundaryField()))
PtrList< surfaceScalarField > phiFs(fluid.phiFs(rAUs))
PtrList< surfaceScalarField > alphafs(phases.size())
PtrList< surfaceScalarField > alpharAUfs(phases.size())
PtrList< volScalarField > rAUs
Definition: pEqn.H:4
volScalarField & alpha
CEqn solve()
multiphaseSystem::phaseModelList & phases
Definition: createFields.H:12
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333