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  forAll(fluid.movingPhases(), movingPhasei)
81  {
82  phaseModel& phase = fluid.movingPhases()[movingPhasei];
83  MRF.correctBoundaryFlux(phase.U(), phase.phiRef());
84  }
85 
86  // Combined buoyancy and force fluxes
87  PtrList<surfaceScalarField> phigFs(phases.size());
88  {
89  const surfaceScalarField ghSnGradRho
90  (
91  "ghSnGradRho",
92  ghf*fvc::snGrad(rho)*mesh.magSf()
93  );
94 
96  {
97  phaseModel& phase = phases[phasei];
98 
99  phigFs.set
100  (
101  phasei,
102  (
104  *(
105  ghSnGradRho
106  - (fvc::interpolate(phase.rho() - rho))*(g & mesh.Sf())
107  - fluid.surfaceTension(phase)*mesh.magSf()
108  )
109  ).ptr()
110  );
111 
112  if (phiFfs.set(phasei))
113  {
114  phigFs[phasei] += phiFfs[phasei];
115  }
116  }
117  }
118 
119  // Predicted fluxes for each phase
120  PtrList<surfaceScalarField> phiHbyAs(phases.size());
121  forAll(fluid.movingPhases(), movingPhasei)
122  {
123  phaseModel& phase = fluid.movingPhases()[movingPhasei];
124 
125  phiHbyAs.set
126  (
127  phase.index(),
129  (
130  IOobject::groupName("phiHbyA", phase.name()),
131  rAUfs[phase.index()]
132  *(
133  fvc::flux(UEqns[phase.index()].H())
134  + alphaRho0fs[phase.index()]
135  *byDt(MRF.absolute(phase.phi()().oldTime()))
136  )
137  - phigFs[phase.index()]
138  )
139  );
140  }
141  fluid.fillFields("phiHbyA", dimForce/dimDensity/dimVelocity, phiHbyAs);
142 
143  // Add explicit drag forces and fluxes
144  PtrList<surfaceScalarField> phiKdPhifs(fluid.phiKdPhifs(rAUfs));
145 
147  {
148  if (phiKdPhifs.set(phasei))
149  {
150  phiHbyAs[phasei] -= phiKdPhifs[phasei];
151  }
152  }
153 
154  // Total predicted flux
156  (
157  IOobject
158  (
159  "phiHbyA",
160  runTime.timeName(),
161  mesh,
162  IOobject::NO_READ,
163  IOobject::AUTO_WRITE
164  ),
165  mesh,
167  );
168 
170  {
171  phiHbyA += alphafs[phasei]*phiHbyAs[phasei];
172  }
173 
174  MRF.makeRelative(phiHbyA);
175 
176  // Construct pressure "diffusivity"
178  (
179  IOobject
180  (
181  "rAUf",
182  runTime.timeName(),
183  mesh
184  ),
185  mesh,
186  dimensionedScalar(dimensionSet(-1, 3, 1, 0, 0), 0)
187  );
188 
190  {
191  rAUf += alphafs[phasei]*alpharAUfs[phasei];
192  }
193 
194  rAUf = mag(rAUf);
195 
196  // Update the fixedFluxPressure BCs to ensure flux consistency
197  {
198  surfaceScalarField::Boundary phib(phi.boundaryField());
199  phib = 0;
201  {
202  phaseModel& phase = phases[phasei];
203  phib +=
204  alphafs[phasei].boundaryField()*phase.phi()().boundaryField();
205  }
206 
207  setSnGrad<fixedFluxPressureFvPatchScalarField>
208  (
209  p_rgh.boundaryFieldRef(),
210  (
211  phiHbyA.boundaryField() - phib
212  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
213  );
214  }
215 
216  // Compressible pressure equations
217  PtrList<fvScalarMatrix> pEqnComps(phases.size());
218  PtrList<volScalarField> dmdts(fluid.dmdts());
220  {
221  phaseModel& phase = phases[phasei];
222  const volScalarField& alpha = phase;
223  volScalarField& rho = phase.thermoRef().rho();
224 
225  pEqnComps.set(phasei, new fvScalarMatrix(p_rgh, dimVolume/dimTime));
226  fvScalarMatrix& pEqnComp = pEqnComps[phasei];
227 
228  // Density variation
229  if (!phase.isochoric() || !phase.pure())
230  {
231  pEqnComp +=
232  (
233  fvc::ddt(alpha, rho) + fvc::div(phase.alphaRhoPhi())
234  - fvc::Sp(fvc::ddt(alpha) + fvc::div(phase.alphaPhi()), rho)
235  )/rho;
236  }
237 
238  // Compressibility
239  if (!phase.incompressible())
240  {
241  if (pimple.transonic())
242  {
243  const surfaceScalarField phid
244  (
245  IOobject::groupName("phid", phase.name()),
246  fvc::interpolate(phase.thermo().psi())*phase.phi()
247  );
248 
249  pEqnComp +=
250  correction
251  (
252  (alpha/rho)*
253  (
254  phase.thermo().psi()*fvm::ddt(p_rgh)
255  + fvm::div(phid, p_rgh)
256  - fvm::Sp(fvc::div(phid), p_rgh)
257  )
258  );
259 
261  (
262  pEqnComps[phasei].faceFluxCorrectionPtr()
263  );
264 
265  pEqnComps[phasei].relax();
266  }
267  else
268  {
269  pEqnComp +=
270  (alpha*phase.thermo().psi()/rho)
271  *correction(fvm::ddt(p_rgh));
272  }
273  }
274 
275  // Option sources
276  if (fvOptions.appliesToField(rho.name()))
277  {
278  pEqnComp -= (fvOptions(alpha, rho) & rho)/rho;
279  }
280 
281  // Mass transfer
282  if (dmdts.set(phasei))
283  {
284  pEqnComp -= dmdts[phasei]/rho;
285  }
286  }
287 
288  // Cache p prior to solve for density update
289  volScalarField p_rgh_0(p_rgh);
290 
291  // Iterate over the pressure equation to correct for non-orthogonality
292  while (pimple.correctNonOrthogonal())
293  {
294  // Construct the transport part of the pressure equation
295  fvScalarMatrix pEqnIncomp
296  (
298  - fvm::laplacian(rAUf, p_rgh)
299  );
300 
301  // Solve
302  {
303  fvScalarMatrix pEqn(pEqnIncomp);
304 
306  {
307  pEqn += pEqnComps[phasei];
308  }
309 
310  if (fluid.incompressible())
311  {
312  pEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
313  }
314 
315  pEqn.solve();
316  }
317 
318  // Correct fluxes and velocities on last non-orthogonal iteration
319  if (pimple.finalNonOrthogonalIter())
320  {
321  phi = phiHbyA + pEqnIncomp.flux();
322 
323  surfaceScalarField mSfGradp("mSfGradp", pEqnIncomp.flux()/rAUf);
324 
325  forAll(fluid.movingPhases(), movingPhasei)
326  {
327  phaseModel& phase = fluid.movingPhases()[movingPhasei];
328 
329  phase.phiRef() =
330  phiHbyAs[phase.index()]
331  + alpharAUfs[phase.index()]*mSfGradp;
332 
333  // Set the phase dilatation rate
334  phase.divU(-pEqnComps[phase.index()] & p_rgh);
335  }
336 
337  if (partialElimination)
338  {
339  fluid.partialEliminationf(rAUfs, alphafs, phiKdPhifs);
340 
341  forAll(fluid.movingPhases(), movingPhasei)
342  {
343  phaseModel& phase = fluid.movingPhases()[movingPhasei];
344 
345  phase.URef() = fvc::reconstruct(MRF.absolute(phase.phi()));
346  phase.URef().correctBoundaryConditions();
347  fvOptions.correct(phase.URef());
348  }
349  }
350  else
351  {
352  forAll(fluid.movingPhases(), movingPhasei)
353  {
354  phaseModel& phase = fluid.movingPhases()[movingPhasei];
355 
356  phase.URef() = fvc::reconstruct(phase.phi());
357  MRF.makeRelative(phase.phiRef());
358  phase.URef().correctBoundaryConditions();
359  fvOptions.correct(phase.URef());
360  }
361  }
362  }
363  }
364 
365  // Update and limit the static pressure
366  p = max(p_rgh + rho*gh, pMin);
367 
368  // Account for static pressure reference
369  if (p_rgh.needReference() && fluid.incompressible())
370  {
372  (
373  "p",
374  p.dimensions(),
376  );
377  }
378 
379  // Limit p_rgh
380  p_rgh = p - rho*gh;
381 
382  // Update densities from change in p_rgh
384  {
385  phaseModel& phase = phases[phasei];
386  phase.thermoRef().rho() += phase.thermo().psi()*(p_rgh - p_rgh_0);
387  }
388 
389  // Correct p_rgh for consistency with p and the updated densities
390  rho = fluid.rho();
391  p_rgh = p - rho*gh;
392  p_rgh.correctBoundaryConditions();
393 }
394 
395 if (!fluid.implicitPhasePressure())
396 {
397  rAUfs.clear();
398 }
forAll(phases, phasei)
Definition: pEqn.H:3
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
rho oldTime()
fv::options & fvOptions
rho
Definition: pEqn.H:1
Info<< "Constructing momentum equations"<< endl;PtrList< fvVectorMatrix > UEqns(phases.size())
pimpleNoLoopControl & pimple
p
Definition: pEqn.H:50
label phasei
Definition: pEqn.H:27
IOMRFZoneList & MRF
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
PtrList< surfaceScalarField > alphafs(phases.size())
volScalarField p_rgh_0(p_rgh)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
p_rgh
Definition: pEqn.H:140
engineTime & runTime
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
phi
Definition: pEqn.H:104
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
phaseSystem & fluid
Definition: createFields.H:11
PtrList< surfaceScalarField > alpharAUfs(phases.size())
dimensionedScalar pMin("pMin", dimPressure, fluid)
const dimensionSet dimForce
const surfaceScalarField & ghf
PtrList< surfaceScalarField > phiFfs(fluid.phiFfs(rAUfs))
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.
const dimensionSet dimDensity
phiHbyA
Definition: pEqn.H:32
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
PtrList< surfaceScalarField > alphaRho0fs(phases.size())
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label pRefCell
Definition: createFields.H:106
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
phaseSystem::phaseModelList & phases
Definition: createFields.H:12
const volScalarField & gh
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
void deleteDemandDrivenData(DataPtr &dataPtr)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
scalar pRefValue
Definition: createFields.H:107
zeroField Sp
Definition: alphaSuSp.H:2
const dimensionSet dimVelocity