pEqn.H
Go to the documentation of this file.
1rho = thermo.rho();
2rho = max(rho, rhoMin);
3rho = min(rho, rhoMax);
4rho.relax();
5
6volScalarField rAU(1.0/UEqn.A());
7surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU));
8volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
9
10if (pimple.nCorrPISO() <= 1)
11{
12 tUEqn.clear();
13}
14
15if (pimple.transonic())
16{
17 surfaceScalarField phid
18 (
19 "phid",
20 fvc::interpolate(psi)
21 *(
22 fvc::flux(HbyA)
23 + MRF.zeroFilter
24 (
25 rhorAUf*fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho)
26 )
27 )
28 );
29
30 MRF.makeRelative(fvc::interpolate(psi), phid);
31
32 while (pimple.correctNonOrthogonal())
33 {
34 fvScalarMatrix pEqn
35 (
36 fvm::ddt(psi, p)
37 + fvm::div(phid, p)
38 - fvm::laplacian(rhorAUf, p)
39 ==
40 parcels.Srho()
41 + fvOptions(psi, p, rho.name())
42 );
43
44 pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
45
46 if (pimple.finalNonOrthogonalIter())
47 {
48 phi == pEqn.flux();
49 }
50 }
51}
52else
53{
54 surfaceScalarField phiHbyA
55 (
56 "phiHbyA",
57 (
58 fvc::flux(rho*HbyA)
59 + rhorAUf*fvc::ddtCorr(rho, U, phi)
60 )
61 );
62
63 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
64
65 // Update the pressure BCs to ensure flux consistency
67
68 while (pimple.correctNonOrthogonal())
69 {
70 fvScalarMatrix pEqn
71 (
72 fvm::ddt(psi, p)
73 + fvc::div(phiHbyA)
74 - fvm::laplacian(rhorAUf, p)
75 ==
76 parcels.Srho()
77 + fvOptions(psi, p, rho.name())
78 );
79
80 pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
81
82 if (pimple.finalNonOrthogonalIter())
83 {
84 phi = phiHbyA + pEqn.flux();
85 }
86 }
87}
88
89#include "rhoEqn.H"
90#include "compressibleContinuityErrs.H"
91
92// Explicitly relax pressure for momentum corrector
93p.relax();
94
95// Recalculate density from the relaxed pressure
96rho = thermo.rho();
97rho = max(rho, rhoMin);
98rho = min(rho, rhoMax);
99rho.relax();
100Info<< "rho min/max : " << min(rho).value() << " " << max(rho).value() << endl;
101
102U = HbyA - rAU*fvc::grad(p);
103U.correctBoundaryConditions();
104fvOptions.correct(U);
105K = 0.5*magSqr(U);
106
107if (thermo.dpdt())
108{
109 dpdt = fvc::ddt(p);
110}
CGAL::Exact_predicates_exact_constructions_kernel K
Y[inertIndex] max(0.0)
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))
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
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
const dimensionedScalar rhoMin
const dimensionedScalar rhoMax
volScalarField & dpdt
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Solve the continuity for density.