pEqn.H
Go to the documentation of this file.
1rho = thermo.rho();
2
3volScalarField rAU(1.0/UEqn.A());
4surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU));
5volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
6
7if (pimple.transonic())
8{
9 surfaceScalarField phid
10 (
11 "phid",
12 fvc::interpolate(psi)
13 *(
14 fvc::flux(HbyA)
15 + MRF.zeroFilter
16 (
17 rhorAUf*fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho)
18 )
19 )
20 );
21
22 MRF.makeRelative(fvc::interpolate(psi), phid);
23
24 while (pimple.correctNonOrthogonal())
25 {
26 fvScalarMatrix pEqn
27 (
28 fvm::ddt(psi, p)
29 + fvm::div(phid, p)
30 - fvm::laplacian(rhorAUf, p)
31 ==
32 coalParcels.Srho()
33 + fvOptions(psi, p, rho.name())
34 );
35
36 pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
37
38 if (pimple.finalNonOrthogonalIter())
39 {
40 phi == pEqn.flux();
41 }
42 }
43}
44else
45{
46 surfaceScalarField phiHbyA
47 (
48 "phiHbyA",
49 (
50 fvc::flux(rho*HbyA)
51 + MRF.zeroFilter(rhorAUf*fvc::ddtCorr(rho, U, phi))
52 )
53 );
54
55 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
56
57 // Update the pressure BCs to ensure flux consistency
59
60 while (pimple.correctNonOrthogonal())
61 {
62 fvScalarMatrix pEqn
63 (
64 fvm::ddt(psi, p)
65 + fvc::div(phiHbyA)
66 - fvm::laplacian(rhorAUf, p)
67 ==
68 coalParcels.Srho()
69 + fvOptions(psi, p, rho.name())
70 );
71
72 pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
73
74 if (pimple.finalNonOrthogonalIter())
75 {
76 phi = phiHbyA + pEqn.flux();
77 }
78 }
79}
80
81#include "rhoEqn.H"
82#include "compressibleContinuityErrs.H"
83
84U = HbyA - rAU*fvc::grad(p);
85U.correctBoundaryConditions();
86fvOptions.correct(U);
87
88K = 0.5*magSqr(U);
89
90if (thermo.dpdt())
91{
92 dpdt = fvc::ddt(p);
93}
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...
U
Definition: pEqn.H:72
volScalarField & p
const volScalarField & psi
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
fvVectorMatrix & UEqn
Definition: UEqn.H:13
phiHbyA
Definition: pcEqn.H:73
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))))
dynamicFvMesh & mesh
volScalarField & dpdt
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
Solve the continuity for density.