pcEqn.H
Go to the documentation of this file.
1if (!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
8const volScalarField psip0(psi*p);
9
10volScalarField rAU(1.0/UEqn.A());
11volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()));
12volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
13
14if (pimple.nCorrPISO() <= 1)
15{
16 tUEqn.clear();
17}
18
19surfaceScalarField phiHbyA
20(
21 "phiHbyA",
22 (
23 fvc::interpolate(rho)*fvc::flux(HbyA)
24 + MRF.zeroFilter
25 (
26 fvc::interpolate(rho*rAU)*fvc::ddtCorr(rho, U, phi, rhoUf)
27 )
28 )
29);
30
31fvc::makeRelative(phiHbyA, rho, U);
32MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
33
34volScalarField rhorAtU("rhorAtU", rho*rAtU);
35
36// Update the pressure BCs to ensure flux consistency
38
39if (pimple.transonic())
40{
41 surfaceScalarField phid
42 (
43 "phid",
44 (fvc::interpolate(psi)/fvc::interpolate(rho))*phiHbyA
45 );
46
47 phiHbyA +=
48 fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf()
49 - fvc::interpolate(psi*p)*phiHbyA/fvc::interpolate(rho);
50
51 HbyA -= (rAU - rAtU)*fvc::grad(p);
52
53 fvScalarMatrix pDDtEqn
54 (
55 fvc::ddt(rho) + psi*correction(fvm::ddt(p))
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}
76else
77{
78 phiHbyA += fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf();
79 HbyA -= (rAU - rAtU)*fvc::grad(p);
80
81 fvScalarMatrix pDDtEqn
82 (
83 fvc::ddt(rho) + psi*correction(fvm::ddt(p))
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#include "rhoEqn.H"
103#include "compressibleContinuityErrs.H"
104
105// Explicitly relax pressure for momentum corrector
106p.relax();
107
108U = HbyA - rAtU*fvc::grad(p);
109U.correctBoundaryConditions();
110fvOptions.correct(U);
111K = 0.5*magSqr(U);
112
114{
115 p.correctBoundaryConditions();
116}
117
118thermo.correctRho(psi*p - psip0, rhoMin, rhoMax) ;
119rho = thermo.rho();
120
121// Correct rhoUf if the mesh is moving
122fvc::correctRhoUf(rhoUf, rho, U, phi);
123
124if (thermo.dpdt())
125{
126 dpdt = fvc::ddt(p);
127
128 if (mesh.moving())
129 {
130 dpdt -= fvc::div(fvc::meshPhi(rho, U), p);
131 }
132}
CGAL::Exact_predicates_exact_constructions_kernel K
fv::options & fvOptions
surfaceScalarField & phi
IOMRFZoneList & MRF
pimpleControl & pimple
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
volScalarField & p
const volScalarField & psi
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
fvVectorMatrix & UEqn
Definition: UEqn.H:13
volScalarField rhorAtU("rhorAtU", rho *rAtU)
phiHbyA
Definition: pcEqn.H:73
U
Definition: pcEqn.H:107
rho
Definition: pcEqn.H:1
volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()))
HbyA
Definition: pcEqn.H:74
const volScalarField psip0(psi *p)
fvScalarMatrix pDDtEqn(fvc::ddt(rho)+psi *correction(fvm::ddt(p))+fvc::div(phiHbyA)==fvOptions(psi, p, rho.name()))
surfaceScalarField phid("phid", fvc::interpolate(psi) *(fvc::flux(HbyA)+MRF.zeroFilter(rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho))))
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > rhoUf
const dimensionedScalar rhoMin
const dimensionedScalar rhoMax
volScalarField & dpdt
const pressureControl & pressureControl
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1