pcEqn.H
Go to the documentation of this file.
1 if (!pimple.SIMPLErho())
2 {
3  rho = thermo.rho();
4 }
5 
6 // Thermodynamic density needs to be updated by psi*d(p) after the
7 // pressure solution
8 const volScalarField psip0(psi*p);
9 
10 volScalarField rAU(1.0/UEqn.A());
11 volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()));
13 
14 if (pimple.nCorrPISO() <= 1)
15 {
16  tUEqn.clear();
17 }
18 
20 (
21  "phiHbyA",
22  (
25  )
26 );
27 
28 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
29 
30 volScalarField rhorAtU("rhorAtU", rho*rAtU);
31 
32 // Update the pressure BCs to ensure flux consistency
33 constrainPressure(p, rho, U, phiHbyA, rhorAtU, MRF);
34 
35 if (pimple.transonic())
36 {
38  (
39  "phid",
41  );
42 
43  phiHbyA +=
44  fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf()
46 
47  HbyA -= (rAU - rAtU)*fvc::grad(p);
48 
50  (
52  + fvc::div(phiHbyA) + fvm::div(phid, p)
53  ==
54  fvOptions(psi, p, rho.name())
55  );
56 
57  while (pimple.correctNonOrthogonal())
58  {
59  fvScalarMatrix pEqn(pDDtEqn - fvm::laplacian(rhorAtU, p));
60 
61  // Relax the pressure equation to ensure diagonal-dominance
62  pEqn.relax();
63 
64  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
65 
66  if (pimple.finalNonOrthogonalIter())
67  {
68  phi = phiHbyA + pEqn.flux();
69  }
70  }
71 }
72 else
73 {
74  phiHbyA += fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf();
75  HbyA -= (rAU - rAtU)*fvc::grad(p);
76 
78  (
80  + fvc::div(phiHbyA)
81  ==
82  fvOptions(psi, p, rho.name())
83  );
84 
85  while (pimple.correctNonOrthogonal())
86  {
87  fvScalarMatrix pEqn(pDDtEqn - fvm::laplacian(rhorAtU, p));
88 
89  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
90 
91  if (pimple.finalNonOrthogonalIter())
92  {
93  phi = phiHbyA + pEqn.flux();
94  }
95  }
96 }
97 
98 // Thermodynamic density update
99 thermo.correctRho(psi*p - psip0);
100 
101 #include "rhoEqn.H"
102 #include "compressibleContinuityErrs.H"
103 
104 // Explicitly relax pressure for momentum corrector
105 p.relax();
106 
107 U = HbyA - rAtU*fvc::grad(p);
108 U.correctBoundaryConditions();
109 fvOptions.correct(U);
110 K = 0.5*magSqr(U);
111 
112 if (pressureControl.limit(p))
113 {
114  p.correctBoundaryConditions();
115  rho = thermo.rho();
116 }
117 else if (pimple.SIMPLErho())
118 {
119  rho = thermo.rho();
120 }
121 
122 if (thermo.dpdt())
123 {
124  dpdt = fvc::ddt(p);
125 }
surfaceScalarField & phi
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
fvScalarMatrix pDDtEqn(fvc::ddt(rho)+psi *correction(fvm::ddt(p))+fvc::div(phiHbyA)==fvOptions(psi, p, rho.name()))
constrainPressure(p, rho, U, phiHbyA, rhorAtU, MRF)
volScalarField rAU(1.0/UEqn.A())
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:170
const dictionary & pimple
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
CGAL::Exact_predicates_exact_constructions_kernel K
HbyA
Definition: pcEqn.H:71
fv::options & fvOptions
phiHbyA
Definition: pcEqn.H:70
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:52
psiReactionThermo & thermo
Definition: createFields.H:31
dynamicFvMesh & mesh
IOMRFZoneList & MRF
volScalarField & dpdt
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
rho
Definition: pcEqn.H:1
dimensioned< scalar > magSqr(const dimensioned< Type > &)
U
Definition: pcEqn.H:104
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.
surfaceScalarField phid("phid", fvc::interpolate(psi) *(fvc::flux(HbyA)+rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho)))
volScalarField rhorAtU("rhorAtU", rho *rAtU)
fvVectorMatrix & UEqn
Definition: UEqn.H:13
const volScalarField & psi
const volScalarField psip0(psi *p)
volScalarField & p
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
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()))