pEqn.H
Go to the documentation of this file.
1 // Face volume fractions
2 PtrList<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
13 rAUs.clear();
14 rAUs.setSize(phases.size());
15 
16 forAll(fluid.movingPhases(), movingPhasei)
17 {
18  phaseModel& phase = fluid.movingPhases()[movingPhasei];
19  const volScalarField& alpha = phase;
20 
21  rAUs.set
22  (
23  phase.index(),
24  new volScalarField
25  (
26  IOobject::groupName("rAU", phase.name()),
27  1.0
28  /(
29  UEqns[phase.index()].A()
30  + byDt(max(phase.residualAlpha() - alpha, scalar(0))*phase.rho())
31  )
32  )
33  );
34 }
35 fluid.fillFields("rAU", dimTime/dimDensity, rAUs);
36 
37 // Phase diagonal coefficients
38 PtrList<surfaceScalarField> alpharAUfs(phases.size());
40 {
41  phaseModel& phase = phases[phasei];
42  const volScalarField& alpha = phase;
43 
44  alpharAUfs.set
45  (
46  phasei,
47  (
48  fvc::interpolate(max(alpha, phase.residualAlpha())*rAUs[phasei])
49  ).ptr()
50  );
51 }
52 
53 // Explicit force fluxes
54 PtrList<surfaceScalarField> phiFs(fluid.phiFs(rAUs));
55 
56 // --- Pressure corrector loop
57 while (pimple.correct())
58 {
59  volScalarField rho("rho", fluid.rho());
60 
61  // Correct p_rgh for consistency with p and the updated densities
62  p_rgh = p - rho*gh;
63 
64  // Correct fixed-flux BCs to be consistent with the velocity BCs
65  fluid.correctBoundaryFlux();
66 
67  // Combined buoyancy and force fluxes
68  PtrList<surfaceScalarField> phigFs(phases.size());
69  {
70  const surfaceScalarField ghSnGradRho
71  (
72  "ghSnGradRho",
73  ghf*fvc::snGrad(rho)*mesh.magSf()
74  );
75 
77  {
78  phaseModel& phase = phases[phasei];
79 
80  phigFs.set
81  (
82  phasei,
83  (
85  *(
86  ghSnGradRho
87  - (fvc::interpolate(phase.rho() - rho))*(g & mesh.Sf())
88  - fluid.surfaceTension(phase)*mesh.magSf()
89  )
90  ).ptr()
91  );
92 
93  if (phiFs.set(phasei))
94  {
95  phigFs[phasei] += phiFs[phasei];
96  }
97  }
98  }
99 
100  // Predicted velocities and fluxes for each phase
101  PtrList<volVectorField> HbyAs(phases.size());
102  PtrList<surfaceScalarField> phiHbyAs(phases.size());
103  {
104  // Correction force fluxes
105  PtrList<surfaceScalarField> ddtCorrByAs(fluid.ddtCorrByAs(rAUs));
106 
107  forAll(fluid.movingPhases(), movingPhasei)
108  {
109  phaseModel& phase = fluid.movingPhases()[movingPhasei];
110  const volScalarField& alpha = phase;
111 
112  HbyAs.set
113  (
114  phase.index(),
116  (
117  rAUs[phase.index()]
118  *(
119  UEqns[phase.index()].H()
120  + byDt
121  (
122  max(phase.residualAlpha() - alpha, scalar(0))
123  *phase.rho()
124  )
125  *phase.U()().oldTime()
126  ),
127  phase.U(),
128  p_rgh
129  )
130  );
131 
132  phiHbyAs.set
133  (
134  phase.index(),
136  (
137  IOobject::groupName("phiHbyA", phase.name()),
138  fvc::flux(HbyAs[phase.index()])
139  - phigFs[phase.index()]
140  - ddtCorrByAs[phase.index()]
141  )
142  );
143  }
144  }
145  fluid.fillFields("HbyA", dimVelocity, HbyAs);
146  fluid.fillFields("phiHbyA", dimForce/dimDensity/dimVelocity, phiHbyAs);
147 
148  // Add explicit drag forces and fluxes
149  PtrList<volVectorField> KdUByAs(fluid.KdUByAs(rAUs));
150  PtrList<surfaceScalarField> phiKdPhis(fluid.phiKdPhis(rAUs));
151 
153  {
154  if (KdUByAs.set(phasei))
155  {
156  HbyAs[phasei] -= KdUByAs[phasei];
157  }
158 
159  if (phiKdPhis.set(phasei))
160  {
161  phiHbyAs[phasei] -= phiKdPhis[phasei];
162  }
163  }
164 
165  // Total predicted flux
167  (
168  IOobject
169  (
170  "phiHbyA",
171  runTime.timeName(),
172  mesh,
173  IOobject::NO_READ,
174  IOobject::AUTO_WRITE
175  ),
176  mesh,
178  );
179 
181  {
182  phiHbyA += alphafs[phasei]*phiHbyAs[phasei];
183  }
184 
185  MRF.makeRelative(phiHbyA);
186  fvc::makeRelative(phiHbyA, fluid.movingPhases()[0].U());
187 
188  // Construct pressure "diffusivity"
190  (
191  IOobject
192  (
193  "rAUf",
194  runTime.timeName(),
195  mesh
196  ),
197  mesh,
198  dimensionedScalar(dimensionSet(-1, 3, 1, 0, 0), 0)
199  );
200 
202  {
203  rAUf += alphafs[phasei]*alpharAUfs[phasei];
204  }
205 
206  rAUf = mag(rAUf);
207 
208  // Update the fixedFluxPressure BCs to ensure flux consistency
209  {
210  surfaceScalarField::Boundary phib
211  (
212  surfaceScalarField::Internal::null(),
213  phi.boundaryField()
214  );
215  phib = 0;
216 
218  {
219  phaseModel& phase = phases[phasei];
220  phib +=
221  alphafs[phasei].boundaryField()*phase.phi()().boundaryField();
222  }
223 
224  setSnGrad<fixedFluxPressureFvPatchScalarField>
225  (
226  p_rgh.boundaryFieldRef(),
227  (
228  phiHbyA.boundaryField() - phib
229  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
230  );
231  }
232 
233  // Compressible pressure equations
234  #include "pEqnComps.H"
235 
236  // Cache p prior to solve for density update
237  volScalarField p_rgh_0(p_rgh);
238 
239  // Iterate over the pressure equation to correct for non-orthogonality
240  while (pimple.correctNonOrthogonal())
241  {
242  // Construct the transport part of the pressure equation
243  fvScalarMatrix pEqnIncomp
244  (
246  - fvm::laplacian(rAUf, p_rgh)
247  );
248 
249  // Solve
250  {
251  fvScalarMatrix pEqn(pEqnIncomp);
252 
254  {
255  pEqn += pEqnComps[phasei];
256  }
257 
258  if (fluid.incompressible())
259  {
260  pEqn.setReference
261  (
262  pressureReference.refCell(),
263  pressureReference.refValue()
264  );
265  }
266 
267  pEqn.solve();
268  }
269 
270  // Correct fluxes and velocities on last non-orthogonal iteration
271  if (pimple.finalNonOrthogonalIter())
272  {
273  phi = phiHbyA + pEqnIncomp.flux();
274 
275  surfaceScalarField mSfGradp("mSfGradp", pEqnIncomp.flux()/rAUf);
276 
277  forAll(fluid.movingPhases(), movingPhasei)
278  {
279  phaseModel& phase = fluid.movingPhases()[movingPhasei];
280 
281  phase.phiRef() =
282  phiHbyAs[phase.index()]
283  + alpharAUfs[phase.index()]*mSfGradp;
284 
285  // Set the phase dilatation rate
286  phase.divU(-pEqnComps[phase.index()] & p_rgh);
287  }
288 
289  // Optionally relax pressure for velocity correction
290  p_rgh.relax();
291 
292  mSfGradp = pEqnIncomp.flux()/rAUf;
293 
294  forAll(fluid.movingPhases(), movingPhasei)
295  {
296  phaseModel& phase = fluid.movingPhases()[movingPhasei];
297 
298  phase.URef() =
299  HbyAs[phase.index()]
301  (
302  alpharAUfs[phase.index()]*mSfGradp
303  - phigFs[phase.index()]
304  );
305  }
306 
307  if (partialElimination)
308  {
309  fluid.partialElimination(rAUs, KdUByAs, alphafs, phiKdPhis);
310  }
311  else
312  {
313  forAll(fluid.movingPhases(), movingPhasei)
314  {
315  phaseModel& phase = fluid.movingPhases()[movingPhasei];
316 
317  MRF.makeRelative(phase.phiRef());
318  fvc::makeRelative(phase.phiRef(), phase.U());
319  }
320  }
321 
322  forAll(fluid.movingPhases(), movingPhasei)
323  {
324  phaseModel& phase = fluid.movingPhases()[movingPhasei];
325 
326  phase.URef().correctBoundaryConditions();
327  phase.correctUf();
328  fvConstraints.constrain(phase.URef());
329  }
330  }
331  }
332 
333  // Update and limit the static pressure
334  p = p_rgh + rho*gh;
336 
337  // Account for static pressure reference
338  if (p_rgh.needReference() && fluid.incompressible())
339  {
341  (
342  "p",
343  p.dimensions(),
344  pressureReference.refValue()
345  - getRefCellValue(p, pressureReference.refCell())
346  );
347  }
348 
349  // Limit p_rgh
350  p_rgh = p - rho*gh;
351 
352  // Update densities from change in p_rgh
354  {
355  phaseModel& phase = phases[phasei];
356  phase.thermoRef().rho() += phase.thermo().psi()*(p_rgh - p_rgh_0);
357  }
358 
359  // Correct p_rgh for consistency with p and the updated densities
360  rho = fluid.rho();
361  p_rgh = p - rho*gh;
362  p_rgh.correctBoundaryConditions();
363 }
364 
365 if (!fluid.implicitPhasePressure())
366 {
367  rAUs.clear();
368 }
forAll(phases, phasei)
Definition: pEqn.H:3
pressureReference & pressureReference
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
pimpleNoLoopControl & pimple
label phasei
Definition: pEqn.H:27
p_rgh
Definition: pEqn.H:160
mixture MRF().makeRelative(phiHbyA)
PtrList< surfaceScalarField > alphafs(phases.size())
UEqns[i]
Definition: UEqn.H:5
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
MRF makeRelative(fvc::interpolate(rho), phiHbyA)
fvMesh & mesh
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
scalar getRefCellValue(const volScalarField &field, const label refCelli)
Return the current value of field in the reference cell.
Definition: findRefCell.C:136
p
Definition: pEqn.H:125
PtrList< fvScalarMatrix > pEqnComps(phases.size())
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionSet dimTime
PtrList< surfaceScalarField > phiFs(fluid.phiFs(rAUs))
phiHbyA
Definition: pEqn.H:30
Foam::fvConstraints & fvConstraints
volScalarField p_rgh_0(p_rgh)
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);if(fluid.found("pMin")){ IOWarningInFunction(fluid)<< "Pressure limits, pMin and pMax, are now read from "<< pimple.dict().name()<< endl;}pressureReference pressureReference(p, p_rgh, pimple.dict(), fluid.incompressible());if(fluid.incompressible()){ p=p_rgh+fluid.rho() *gh;}if(p_rgh.needReference() &&fluid.incompressible()){ p+=dimensionedScalar("p", p.dimensions(), pressureReference.refValue() - getRefCellValue(p, pressureReference.refCell()));}p_rgh=p - fluid.rho() *gh;mesh.schemes().setFluxRequired(p_rgh.name());PtrList< volScalarField > rAUs
Definition: createFields.H:67
const dimensionSet dimFlux
const dimensionSet dimDensity
phaseSystem & fluid
Definition: createFields.H:11
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
const dimensionSet dimForce
PtrList< surfaceScalarField > alpharAUfs(phases.size())
const surfaceScalarField & ghf
const dimensionSet dimVelocity
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
tmp< volScalarField > byDt(const volScalarField &vf)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
rho oldTime()
phi
Definition: pEqn.H:98
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
phaseSystem::phaseModelList & phases
Definition: createFields.H:12
const volScalarField & gh
dimensioned< scalar > mag(const dimensioned< Type > &)
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
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
const dimensionedVector & g
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
rho
Definition: pEqn.H:1