pEqnComps.H
Go to the documentation of this file.
1 PtrList<fvScalarMatrix> pEqnComps(phases.size());
2 
3 {
4  PtrList<volScalarField> dmdts(fluid.dmdts());
5  PtrList<volScalarField> d2mdtdps(fluid.d2mdtdps());
6 
8  {
9  phaseModel& phase = phases[phasei];
10  const volScalarField& alpha = phase;
11  volScalarField& rho = phase.thermoRef().rho();
12 
14  fvScalarMatrix& pEqnComp = pEqnComps[phasei];
15 
16  // Density variation
17  if (!phase.isochoric() || !phase.pure())
18  {
19  pEqnComp +=
20  (
21  fvc::ddt(alpha, rho) + fvc::div(phase.alphaRhoPhi())
22  - fvc::Sp(fvc::ddt(alpha) + fvc::div(phase.alphaPhi()), rho)
23  )/rho;
24  }
25 
26  // Mesh dilatation correction
27  if (mesh.moving())
28  {
29  pEqnComp += fvc::div(mesh.phi())*alpha;
30  }
31 
32  // Compressibility
33  if (!phase.incompressible())
34  {
35  if (pimple.transonic())
36  {
37  const surfaceScalarField phid
38  (
39  IOobject::groupName("phid", phase.name()),
40  fvc::interpolate(phase.thermo().psi())*phase.phi()
41  );
42 
43  pEqnComp +=
45  (
46  (alpha/rho)*
47  (
48  phase.thermo().psi()*fvm::ddt(p_rgh)
49  + fvm::div(phid, p_rgh)
50  - fvm::Sp(fvc::div(phid), p_rgh)
51  )
52  );
53 
54  pEqnComps[phasei].relax();
55  }
56  else
57  {
58  pEqnComp +=
59  (alpha*phase.thermo().psi()/rho)
61  }
62  }
63 
64  // Option sources
65  if (fvModels.addsSupToField(rho.name()))
66  {
67  pEqnComp -= (fvModels.source(alpha, rho) & rho)/rho;
68  }
69 
70  // Mass transfer
71  if (dmdts.set(phasei))
72  {
73  pEqnComp -= dmdts[phasei]/rho;
74  }
75  if (d2mdtdps.set(phasei))
76  {
77  pEqnComp -= correction(fvm::Sp(d2mdtdps[phasei]/rho, p_rgh));
78  }
79  }
80 }
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
volScalarField & p_rgh
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
pimpleNoLoopControl & pimple
label phasei
Definition: pEqn.H:27
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))
fvMesh & mesh
PtrList< fvScalarMatrix > pEqnComps(phases.size())
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:58
const dimensionSet dimTime
phaseSystem & fluid
Definition: createFields.H:11
forAll(phases, phasei)
Definition: pEqnComps.H:7
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 tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
Foam::fvModels & fvModels
const dimensionSet dimVolume
phaseSystem::phaseModelList & phases
Definition: createFields.H:12
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual bool addsSupToField(const word &fieldName) const
Return true if an fvModel adds a source term to the given.
Definition: fvModels.C:224
PtrList< volScalarField > d2mdtdps(fluid.d2mdtdps())