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  (
24  + MRF.zeroFilter
25  (
27  )
28  )
29 );
30 
31 fvc::makeRelative(phiHbyA, rho, U);
32 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
33 
34 volScalarField rhorAtU("rhorAtU", rho*rAtU);
35 
36 // Update the pressure BCs to ensure flux consistency
37 constrainPressure(p, rho, U, phiHbyA, rhorAtU, MRF);
38 
39 if (pimple.transonic())
40 {
42  (
43  "phid",
45  );
46 
47  phiHbyA +=
48  fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf()
50 
51  HbyA -= (rAU - rAtU)*fvc::grad(p);
52 
54  (
56  + fvc::div(phiHbyA) + fvm::div(phid, p)
57  ==
58  fvOptions(psi, p, rho.name())
59  );
60 
61  while (pimple.correctNonOrthogonal())
62  {
63  fvScalarMatrix pEqn(pDDtEqn - fvm::laplacian(rhorAtU, p));
64 
65  // Relax the pressure equation to ensure diagonal-dominance
66  pEqn.relax();
67 
68  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
69 
70  if (pimple.finalNonOrthogonalIter())
71  {
72  phi = phiHbyA + pEqn.flux();
73  }
74  }
75 }
76 else
77 {
78  phiHbyA += fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf();
79  HbyA -= (rAU - rAtU)*fvc::grad(p);
80 
82  (
84  + fvc::div(phiHbyA)
85  ==
86  fvOptions(psi, p, rho.name())
87  );
88 
89  while (pimple.correctNonOrthogonal())
90  {
91  fvScalarMatrix pEqn(pDDtEqn - fvm::laplacian(rhorAtU, p));
92 
93  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
94 
95  if (pimple.finalNonOrthogonalIter())
96  {
97  phi = phiHbyA + pEqn.flux();
98  }
99  }
100 }
101 
102 // Thermodynamic density update
103 thermo.correctRho(psi*p - psip0);
104 
105 #include "rhoEqn.H"
106 #include "compressibleContinuityErrs.H"
107 
108 // Explicitly relax pressure for momentum corrector
109 p.relax();
110 
111 U = HbyA - rAtU*fvc::grad(p);
112 U.correctBoundaryConditions();
113 fvOptions.correct(U);
114 K = 0.5*magSqr(U);
115 
116 if (pressureControl.limit(p))
117 {
118  p.correctBoundaryConditions();
119  rho = thermo.rho();
120 }
121 else if (pimple.simpleRho())
122 {
123  rho = thermo.rho();
124 }
125 
126 // Correct rhoUf if the mesh is moving
128 
129 if (thermo.dpdt())
130 {
131  dpdt = fvc::ddt(p);
132 
133  if (mesh.moving())
134  {
135  dpdt -= fvc::div(fvc::meshPhi(rho, U), p);
136  }
137 }
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()))
fv::options & fvOptions
pimpleNoLoopControl & pimple
surfaceScalarField & phi
IOMRFZoneList & MRF
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
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
rhoReactionThermo & thermo
Definition: createFields.H:28
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:74
surfaceScalarField phid("phid", fvc::interpolate(psi) *(fvc::flux(HbyA)+MRF.zeroFilter(rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho))))
phiHbyA
Definition: pcEqn.H:73
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
volScalarField & dpdt
dynamicFvMesh & mesh
MRF makeRelative(fvc::interpolate(rho), phiHbyA)
pressureControl & pressureControl
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:107
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.
volScalarField rhorAtU("rhorAtU", rho *rAtU)
fvVectorMatrix & UEqn
Definition: UEqn.H:13
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:34
void correctRhoUf(autoPtr< surfaceVectorField > &rhoUf, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi)
Definition: fvcMeshPhi.C:241
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
autoPtr< surfaceVectorField > rhoUf
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()))