pEqn.H
Go to the documentation of this file.
1 // Face volume fractions
2 PtrList<surfaceScalarField> alphafs(phases.size());
3 PtrList<surfaceScalarField> alphaRho0fs(phases.size());
5 {
6  phaseModel& phase = phases[phasei];
7  const volScalarField& alpha = phase;
8 
10  alphafs[phasei].rename("pEqn" + alphafs[phasei].name());
11 
12  alphaRho0fs.set
13  (
14  phasei,
15  (
17  (
18  max(alpha.oldTime(), phase.residualAlpha())
19  *phase.rho()().oldTime()
20  )
21  ).ptr()
22  );
23 }
24 
25 // Diagonal coefficients
26 PtrList<surfaceScalarField> rAUfs(phases.size());
27 {
28  PtrList<surfaceScalarField> AFfs(fluid.AFfs());
29 
30  forAll(fluid.movingPhases(), movingPhasei)
31  {
32  phaseModel& phase = fluid.movingPhases()[movingPhasei];
33 
34  rAUfs.set
35  (
36  phase.index(),
38  (
39  IOobject::groupName("rAUf", phase.name()),
40  1.0
41  /(
42  byDt(alphaRho0fs[phase.index()])
43  + fvc::interpolate(UEqns[phase.index()].A())
44  + AFfs[phase.index()]
45  )
46  )
47  );
48  }
49 }
50 fluid.fillFields("rAUf", dimTime/dimDensity, rAUfs);
51 
52 // Phase diagonal coefficients
53 PtrList<surfaceScalarField> alpharAUfs(phases.size());
55 {
56  phaseModel& phase = phases[phasei];
57  alpharAUfs.set
58  (
59  phase.index(),
60  (
61  max(alphafs[phase.index()], phase.residualAlpha())
62  *rAUfs[phase.index()]
63  ).ptr()
64  );
65 }
66 
67 // Explicit force fluxes
68 PtrList<surfaceScalarField> phiFfs(fluid.phiFfs(rAUfs));
69 
70 // --- Pressure corrector loop
71 while (pimple.correct())
72 {
73  volScalarField rho("rho", fluid.rho());
74 
75  // Correct p_rgh for consistency with p and the updated densities
76  p_rgh = p - rho*gh;
77 
78  // Correct fixed-flux BCs to be consistent with the velocity BCs
79  forAll(fluid.movingPhases(), movingPhasei)
80  {
81  phaseModel& phase = fluid.movingPhases()[movingPhasei];
82  MRF.correctBoundaryFlux(phase.U(), phase.phiRef());
83  }
84 
85  // Combined buoyancy and force fluxes
86  PtrList<surfaceScalarField> phigFs(phases.size());
87  {
89  (
90  "ghSnGradRho",
91  ghf*fvc::snGrad(rho)*mesh.magSf()
92  );
93 
95  {
96  phaseModel& phase = phases[phasei];
97 
98  phigFs.set
99  (
100  phasei,
101  (
103  *(
105  - (fvc::interpolate(phase.rho() - rho))*(g & mesh.Sf())
106  - fluid.surfaceTension(phase)*mesh.magSf()
107  )
108  ).ptr()
109  );
110 
111  if (phiFfs.set(phasei))
112  {
113  phigFs[phasei] += phiFfs[phasei];
114  }
115  }
116  }
117 
118  // Predicted fluxes for each phase
119  PtrList<surfaceScalarField> phiHbyAs(phases.size());
120  forAll(fluid.movingPhases(), movingPhasei)
121  {
122  phaseModel& phase = fluid.movingPhases()[movingPhasei];
123 
124  phiHbyAs.set
125  (
126  phase.index(),
128  (
129  IOobject::groupName("phiHbyA", phase.name()),
130  rAUfs[phase.index()]
131  *(
132  fvc::flux(UEqns[phase.index()].H())
133  + alphaRho0fs[phase.index()]
134  *byDt(MRF.absolute(phase.phi()().oldTime()))
135  )
136  - phigFs[phase.index()]
137  )
138  );
139  }
140  fluid.fillFields("phiHbyA", dimForce/dimDensity/dimVelocity, phiHbyAs);
141 
142  // Add explicit drag forces and fluxes if not doing partial elimination
143  if (!partialElimination)
144  {
145  PtrList<surfaceScalarField> phiKdPhifs(fluid.phiKdPhifs(rAUfs));
146 
148  {
149  if (phiKdPhifs.set(phasei))
150  {
151  phiHbyAs[phasei] -= phiKdPhifs[phasei];
152  }
153  }
154  }
155 
156  // Total predicted flux
158  (
159  IOobject
160  (
161  "phiHbyA",
162  runTime.timeName(),
163  mesh,
164  IOobject::NO_READ,
165  IOobject::AUTO_WRITE
166  ),
167  mesh,
169  );
170 
172  {
174  }
175 
176  // Add explicit drag fluxes if doing partial elimination
177  if (partialElimination)
178  {
179  PtrList<surfaceScalarField> phiKdPhifs(fluid.phiKdPhifs(rAUfs));
180 
182  {
183  if (phiKdPhifs.set(phasei))
184  {
185  phiHbyA -= alphafs[phasei]*phiKdPhifs[phasei];
186  }
187  }
188  }
189 
190  MRF.makeRelative(phiHbyA);
191 
192  // Construct pressure "diffusivity"
194  (
195  IOobject
196  (
197  "rAUf",
198  runTime.timeName(),
199  mesh
200  ),
201  mesh,
202  dimensionedScalar("rAUf", dimensionSet(-1, 3, 1, 0, 0), 0)
203  );
204 
206  {
208  }
209 
210  rAUf = mag(rAUf);
211 
212  // Update the fixedFluxPressure BCs to ensure flux consistency
213  {
214  surfaceScalarField::Boundary phib(phi.boundaryField());
215  phib = 0;
217  {
218  phaseModel& phase = phases[phasei];
219  phib +=
220  alphafs[phasei].boundaryField()*phase.phi()().boundaryField();
221  }
222 
224  (
225  p_rgh.boundaryFieldRef(),
226  (
227  phiHbyA.boundaryField() - phib
228  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
229  );
230  }
231 
232  // Compressible pressure equations
233  PtrList<fvScalarMatrix> pEqnComps(phases.size());
234  PtrList<volScalarField> dmdts(fluid.dmdts());
236  {
237  phaseModel& phase = phases[phasei];
238  const volScalarField& alpha = phase;
239  volScalarField& rho = phase.thermoRef().rho();
240 
241  if (phase.compressible())
242  {
243  if (pimple.transonic())
244  {
246  (
247  IOobject::groupName("phid", phase.name()),
248  fvc::interpolate(phase.thermo().psi())*phase.phi()
249  );
250 
251  pEqnComps.set
252  (
253  phasei,
254  (
255  (
256  fvc::ddt(alpha, rho) + fvc::div(phase.alphaRhoPhi())
257  - fvc::Sp
258  (
259  fvc::ddt(alpha) + fvc::div(phase.alphaPhi()),
260  rho
261  )
262  )/rho
263  + correction
264  (
265  (alpha/rho)*
266  (
267  phase.thermo().psi()*fvm::ddt(p_rgh)
268  + fvm::div(phid, p_rgh)
270  )
271  )
272  ).ptr()
273  );
274 
275  pEqnComps[phasei].faceFluxCorrectionPtr(nullptr);
276 
277  pEqnComps[phasei].relax();
278  }
279  else
280  {
281  pEqnComps.set
282  (
283  phasei,
284  (
285  (
286  fvc::ddt(alpha, rho) + fvc::div(phase.alphaRhoPhi())
287  - fvc::Sp
288  (
289  (fvc::ddt(alpha) + fvc::div(phase.alphaPhi())),
290  rho
291  )
292  )/rho
293  + (alpha*phase.thermo().psi()/rho)
295  ).ptr()
296  );
297  }
298  }
299 
300  if (fvOptions.appliesToField(rho.name()))
301  {
302  tmp<fvScalarMatrix> optEqn = fvOptions(alpha, rho);
303  if (pEqnComps.set(phasei))
304  {
305  pEqnComps[phasei] -= (optEqn&rho)/rho;
306  }
307  else
308  {
309  pEqnComps.set
310  (
311  phasei,
312  fvm::Su(- (optEqn&rho)/rho, p_rgh).ptr()
313  );
314  }
315  }
316 
317  if (dmdts.set(phasei))
318  {
319  if (pEqnComps.set(phasei))
320  {
321  pEqnComps[phasei] -= dmdts[phasei]/rho;
322  }
323  else
324  {
325  pEqnComps.set
326  (
327  phasei,
328  fvm::Su(- dmdts[phasei]/rho, p_rgh)
329  );
330  }
331  }
332  }
333 
334  // Cache p prior to solve for density update
336 
337  // Iterate over the pressure equation to correct for non-orthogonality
338  while (pimple.correctNonOrthogonal())
339  {
340  // Construct the transport part of the pressure equation
341  fvScalarMatrix pEqnIncomp
342  (
345  );
346 
347  {
348  fvScalarMatrix pEqn(pEqnIncomp);
349 
351  {
352  if (pEqnComps.set(phasei))
353  {
354  pEqn += pEqnComps[phasei];
355  }
356  }
357 
358  pEqn.solve(p_rgh.select(pimple.finalInnerIter()));
359  }
360 
361  // Correct fluxes and velocities on last non-orthogonal iteration
362  if (pimple.finalNonOrthogonalIter())
363  {
364  phi = phiHbyA + pEqnIncomp.flux();
365 
366  surfaceScalarField mSfGradp("mSfGradp", pEqnIncomp.flux()/rAUf);
367 
368  forAll(fluid.movingPhases(), movingPhasei)
369  {
370  phaseModel& phase = fluid.movingPhases()[movingPhasei];
371 
372  phase.phiRef() =
373  phiHbyAs[phase.index()]
374  + alpharAUfs[phase.index()]*mSfGradp;
375 
376  // Set the phase dilatation rates
377  if (pEqnComps.set(phase.index()))
378  {
379  phase.divU(-pEqnComps[phase.index()] & p_rgh);
380  }
381  }
382 
383  if (partialElimination)
384  {
385  fluid.partialEliminationf(rAUfs);
386  }
387 
388  // Optionally relax pressure for velocity correction
389  p_rgh.relax();
390 
391  mSfGradp = pEqnIncomp.flux()/rAUf;
392 
393  forAll(fluid.movingPhases(), movingPhasei)
394  {
395  phaseModel& phase = fluid.movingPhases()[movingPhasei];
396 
397  phase.URef() = fvc::reconstruct(MRF.absolute(phase.phi()));
398  phase.URef().correctBoundaryConditions();
399  fvOptions.correct(phase.URef());
400  }
401  }
402  }
403 
404  // Update and limit the static pressure
405  p = max(p_rgh + rho*gh, pMin);
406 
407  // Limit p_rgh
408  p_rgh = p - rho*gh;
409 
410  // Update densities from change in p_rgh
412  {
413  phaseModel& phase = phases[phasei];
414  phase.thermoRef().rho() += phase.thermo().psi()*(p_rgh - p_rgh_0);
415  }
416 
417  // Correct p_rgh for consistency with p and the updated densities
418  rho = fluid.rho();
419  p_rgh = p - rho*gh;
420  p_rgh.correctBoundaryConditions();
421 }
PtrList< surfaceScalarField > AFfs(fluid.AFfs())
twoPhaseSystem & fluid
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
multiphaseSystem::phaseModelList & phases
Definition: createFields.H:12
Info<< "Creating field kinetic energy K\"<< 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()
Definition: createK.H:12
rho
Definition: pEqn.H:1
zeroField Su
Definition: alphaSuSp.H:1
PtrList< surfaceScalarField > alphaRho0fs(phases.size())
p
Definition: pEqn.H:50
label phasei
Definition: pEqn.H:27
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionedScalar & pMin
phiHbyA
Definition: pEqn.H:20
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
forAll(U.boundaryField(), patchi)
Definition: pEqn.H:18
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
engineTime & runTime
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:40
PtrList< surfaceScalarField > alpharAUfs(phases.size())
surfaceScalarField phid("phid", fvc::interpolate(psi) *(fvc::flux(HbyA)+MRF.zeroFilter(rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho))))
fv::options & fvOptions
const surfaceScalarField & ghf
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
IOMRFZoneList & MRF
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
PtrList< surfaceScalarField > rAUfs
Definition: pEqn.H:30
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:436
p_rgh
Definition: pEqn.H:139
const uniformDimensionedVectorField & g
const dimensionSet dimForce
pimpleControl & pimple
const dimensionSet dimDensity
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
surfaceScalarField ghSnGradRho(ghf *fvc::snGrad(rho) *mesh.magSf())
PtrList< fvVectorMatrix > UEqns(fluid.phases().size())
PtrList< surfaceScalarField > alphafs(phases.size())
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
mesh interpolate(rAU)
setSnGrad< fixedFluxPressureFvPatchScalarField >(p_rgh.boundaryFieldRef(),(phiHbyA.boundaryField() - MRF.relative(phib))/(mesh.magSf().boundaryField() *rAUf.boundaryField()))
const volScalarField & gh
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
volScalarField p_rgh_0(p_rgh)
PtrList< surfaceScalarField > phiHbyAs(fluid.phases().size())
PtrList< surfaceScalarField > phiFfs(fluid.phiFfs(rAUfs))
phi
Definition: pEqn.H:18
phib
Definition: pEqn.H:190
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
zeroField Sp
Definition: alphaSuSp.H:2
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity