pEqn.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());
13 
14 if (pimple.nCorrPISO() <= 1)
15 {
16  tUEqn.clear();
17 }
18 
20 (
21  "phiHbyA",
23  + MRF.zeroFilter(rhorAUf*fvc::ddtCorr(rho, U, phi, rhoUf))
24 );
25 
27 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
28 
29 // Update the pressure BCs to ensure flux consistency
30 constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF);
31 
32 if (pimple.transonic())
33 {
35  (
36  "phid",
38  );
39 
40  phiHbyA -= fvc::interpolate(psi*p)*phiHbyA/fvc::interpolate(rho);
41 
43  (
45  + fvc::div(phiHbyA) + fvm::div(phid, p)
46  ==
47  fvOptions(psi, p, rho.name())
48  );
49 
50  while (pimple.correctNonOrthogonal())
51  {
52  fvScalarMatrix pEqn(pDDtEqn - fvm::laplacian(rhorAUf, p));
53 
54  // Relax the pressure equation to ensure diagonal-dominance
55  pEqn.relax();
56 
57  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
58 
59  if (pimple.finalNonOrthogonalIter())
60  {
61  phi = phiHbyA + pEqn.flux();
62  }
63  }
64 }
65 else
66 {
68  (
70  + fvc::div(phiHbyA)
71  ==
72  fvOptions(psi, p, rho.name())
73  );
74 
75  while (pimple.correctNonOrthogonal())
76  {
77  fvScalarMatrix pEqn(pDDtEqn - fvm::laplacian(rhorAUf, p));
78 
79  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
80 
81  if (pimple.finalNonOrthogonalIter())
82  {
83  phi = phiHbyA + pEqn.flux();
84  }
85  }
86 }
87 
88 // Thermodynamic density update
89 thermo.correctRho(psi*p - psip0);
90 
91 #include "rhoEqn.H"
92 #include "compressibleContinuityErrs.H"
93 
94 // Explicitly relax pressure for momentum corrector
95 p.relax();
96 
97 U = HbyA - rAU*fvc::grad(p);
98 U.correctBoundaryConditions();
99 fvOptions.correct(U);
100 K = 0.5*magSqr(U);
101 
102 if (pressureControl.limit(p))
103 {
104  p.correctBoundaryConditions();
105  rho = thermo.rho();
106 }
107 else if (pimple.simpleRho())
108 {
109  rho = thermo.rho();
110 }
111 
112 // Correct rhoUf if the mesh is moving
114 
115 if (thermo.dpdt())
116 {
117  dpdt = fvc::ddt(p);
118 
119  if (mesh.moving())
120  {
121  dpdt -= fvc::div(fvc::meshPhi(rho, U), p);
122  }
123 }
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
rho
Definition: pEqn.H:1
pimpleNoLoopControl & pimple
p
Definition: pEqn.H:50
volScalarField rAU(1.0/UEqn.A())
IOMRFZoneList & MRF
rhoUf
Definition: pEqn.H:94
phiHbyA
Definition: pEqn.H:20
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
surfaceScalarField phid("phid", fvc::interpolate(psi) *(fvc::flux(HbyA)+MRF.zeroFilter(rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho))))
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
pressureControl & pressureControl
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
dimensioned< scalar > magSqr(const dimensioned< Type > &)
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.
U
Definition: pEqn.H:72
volVectorField & HbyA
Definition: pEqn.H:13
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
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
phi
Definition: pEqn.H:18
const volScalarField & psi
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
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
MRF makeRelative(fvc::interpolate(rho), phiHbyA)
const volScalarField psip0(psi *p)