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 
9  alphafs.set(phasei, fvc::interpolate(alpha).ptr());
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 rAUfs.clear();
27 rAUfs.setSize(phases.size());
28 {
29  PtrList<surfaceScalarField> AFfs(fluid.AFfs());
30 
31  forAll(fluid.movingPhases(), movingPhasei)
32  {
33  phaseModel& phase = fluid.movingPhases()[movingPhasei];
34 
35  rAUfs.set
36  (
37  phase.index(),
39  (
40  IOobject::groupName("rAUf", phase.name()),
41  1.0
42  /(
43  byDt(alphaRho0fs[phase.index()])
44  + fvc::interpolate(UEqns[phase.index()].A())
45  + AFfs[phase.index()]
46  )
47  )
48  );
49  }
50 }
51 fluid.fillFields("rAUf", dimTime/dimDensity, rAUfs);
52 
53 // Phase diagonal coefficients
54 PtrList<surfaceScalarField> alpharAUfs(phases.size());
56 {
57  phaseModel& phase = phases[phasei];
58  alpharAUfs.set
59  (
60  phase.index(),
61  (
62  max(alphafs[phase.index()], phase.residualAlpha())
63  *rAUfs[phase.index()]
64  ).ptr()
65  );
66 }
67 
68 // Explicit force fluxes
69 PtrList<surfaceScalarField> phiFfs(fluid.phiFfs(rAUfs));
70 
71 // --- Pressure corrector loop
72 while (pimple.correct())
73 {
74  volScalarField rho("rho", fluid.rho());
75 
76  // Correct p_rgh for consistency with p and the updated densities
77  p_rgh = p - rho*gh;
78 
79  // Correct fixed-flux BCs to be consistent with the velocity BCs
80  fluid.correctBoundaryFlux();
81 
82  // Combined buoyancy and force fluxes
83  PtrList<surfaceScalarField> phigFs(phases.size());
84  {
85  const surfaceScalarField ghSnGradRho
86  (
87  "ghSnGradRho",
88  ghf*fvc::snGrad(rho)*mesh.magSf()
89  );
90 
92  {
93  phaseModel& phase = phases[phasei];
94 
95  phigFs.set
96  (
97  phasei,
98  (
100  *(
101  ghSnGradRho
102  - (fvc::interpolate(phase.rho() - rho))*(g & mesh.Sf())
103  - fluid.surfaceTension(phase)*mesh.magSf()
104  )
105  ).ptr()
106  );
107 
108  if (phiFfs.set(phasei))
109  {
110  phigFs[phasei] += phiFfs[phasei];
111  }
112  }
113  }
114 
115  // Predicted fluxes for each phase
116  PtrList<surfaceScalarField> phiHbyAs(phases.size());
117  forAll(fluid.movingPhases(), movingPhasei)
118  {
119  phaseModel& phase = fluid.movingPhases()[movingPhasei];
120 
121  phiHbyAs.set
122  (
123  phase.index(),
125  (
126  rAUfs[phase.index()]
127  *(
128  fvc::flux(UEqns[phase.index()].H())
129  + alphaRho0fs[phase.index()]
130  *byDt
131  (
132  phase.Uf().valid()
133  ? (mesh.Sf() & phase.Uf()().oldTime())
134  : MRF.absolute(phase.phi()().oldTime())
135  )
136  )
137  - phigFs[phase.index()],
138  phase.U(),
139  p_rgh
140  )
141  );
142  }
143  fluid.fillFields("phiHbyA", dimForce/dimDensity/dimVelocity, phiHbyAs);
144 
145  // Add explicit drag forces and fluxes
146  PtrList<surfaceScalarField> phiKdPhifs(fluid.phiKdPhifs(rAUfs));
147 
149  {
150  if (phiKdPhifs.set(phasei))
151  {
152  phiHbyAs[phasei] -= phiKdPhifs[phasei];
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  {
173  phiHbyA += alphafs[phasei]*phiHbyAs[phasei];
174  }
175 
176  MRF.makeRelative(phiHbyA);
178 
179  // Construct pressure "diffusivity"
181  (
182  IOobject
183  (
184  "rAUf",
185  runTime.timeName(),
186  mesh
187  ),
188  mesh,
189  dimensionedScalar(dimensionSet(-1, 3, 1, 0, 0), 0)
190  );
191 
193  {
194  rAUf += alphafs[phasei]*alpharAUfs[phasei];
195  }
196 
197  rAUf = mag(rAUf);
198 
199  // Update the fixedFluxPressure BCs to ensure flux consistency
200  {
201  surfaceScalarField::Boundary phib
202  (
203  surfaceScalarField::Internal::null(),
204  phi.boundaryField()
205  );
206  phib = 0;
207 
209  {
210  phaseModel& phase = phases[phasei];
211  phib +=
212  alphafs[phasei].boundaryField()*phase.phi()().boundaryField();
213  }
214 
215  setSnGrad<fixedFluxPressureFvPatchScalarField>
216  (
217  p_rgh.boundaryFieldRef(),
218  (
219  phiHbyA.boundaryField() - phib
220  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
221  );
222  }
223 
224  // Compressible pressure equations
225  #include "pEqnComps.H"
226 
227  // Cache p prior to solve for density update
228  volScalarField p_rgh_0(p_rgh);
229 
230  // Iterate over the pressure equation to correct for non-orthogonality
231  while (pimple.correctNonOrthogonal())
232  {
233  // Construct the transport part of the pressure equation
234  fvScalarMatrix pEqnIncomp
235  (
237  - fvm::laplacian(rAUf, p_rgh)
238  );
239 
240  // Solve
241  {
242  fvScalarMatrix pEqn(pEqnIncomp);
243 
245  {
246  pEqn += pEqnComps[phasei];
247  }
248 
249  if (fluid.incompressible())
250  {
251  pEqn.setReference
252  (
253  pressureReference.refCell(),
254  pressureReference.refValue()
255  );
256  }
257 
258  pEqn.solve();
259  }
260 
261  // Correct fluxes and velocities on last non-orthogonal iteration
262  if (pimple.finalNonOrthogonalIter())
263  {
264  phi = phiHbyA + pEqnIncomp.flux();
265 
266  surfaceScalarField mSfGradp("mSfGradp", pEqnIncomp.flux()/rAUf);
267 
268  forAll(fluid.movingPhases(), movingPhasei)
269  {
270  phaseModel& phase = fluid.movingPhases()[movingPhasei];
271 
272  phase.phiRef() =
273  phiHbyAs[phase.index()]
274  + alpharAUfs[phase.index()]*mSfGradp;
275 
276  // Set the phase dilatation rate
277  phase.divU(-pEqnComps[phase.index()] & p_rgh);
278  }
279 
280  if (partialElimination)
281  {
282  fluid.partialEliminationf(rAUfs, alphafs, phiKdPhifs);
283  }
284  else
285  {
286  forAll(fluid.movingPhases(), movingPhasei)
287  {
288  phaseModel& phase = fluid.movingPhases()[movingPhasei];
289 
290  MRF.makeRelative(phase.phiRef());
291  fvc::makeRelative(phase.phiRef(), phase.U());
292  }
293  }
294 
295  forAll(fluid.movingPhases(), movingPhasei)
296  {
297  phaseModel& phase = fluid.movingPhases()[movingPhasei];
298 
299  phase.URef() = fvc::reconstruct
300  (
301  fvc::absolute(MRF.absolute(phase.phi()), phase.U())
302  );
303 
304  phase.URef().correctBoundaryConditions();
305  phase.correctUf();
306  fvConstraints.constrain(phase.URef());
307  }
308  }
309  }
310 
311  // Update and limit the static pressure
312  p = p_rgh + rho*gh;
314 
315  // Account for static pressure reference
316  if (p_rgh.needReference() && fluid.incompressible())
317  {
319  (
320  "p",
321  p.dimensions(),
322  pressureReference.refValue()
323  - getRefCellValue(p, pressureReference.refCell())
324  );
325  }
326 
327  // Limit p_rgh
328  p_rgh = p - rho*gh;
329 
330  // Update densities from change in p_rgh
332  {
333  phaseModel& phase = phases[phasei];
334  phase.thermoRef().rho() += phase.thermo().psi()*(p_rgh - p_rgh_0);
335  }
336 
337  // Correct p_rgh for consistency with p and the updated densities
338  rho = fluid.rho();
339  p_rgh = p - rho*gh;
340  p_rgh.correctBoundaryConditions();
341 }
342 
343 if (!fluid.implicitPhasePressure())
344 {
345  rAUfs.clear();
346 }
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
U
Definition: pEqn.H:72
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)
tmp< surfaceScalarField > constrainPhiHbyA(const tmp< surfaceScalarField > &tphiHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:79
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
phiHbyA
Definition: pEqn.H:30
Foam::fvConstraints & fvConstraints
volScalarField p_rgh_0(p_rgh)
const dimensionSet dimFlux
const dimensionSet dimDensity
phaseSystem & fluid
Definition: createFields.H:11
const dimensionSet dimForce
PtrList< surfaceScalarField > alpharAUfs(phases.size())
const surfaceScalarField & ghf
PtrList< surfaceScalarField > phiFfs(fluid.phiFfs(rAUfs))
const dimensionSet dimVelocity
PtrList< surfaceScalarField > rAUfs
Definition: createFields.H:68
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)
PtrList< surfaceScalarField > alphaRho0fs(phases.size())
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
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