pEqn.H
Go to the documentation of this file.
1dimensionedScalar compressibility = fvc::domainIntegrate(psi);
2bool compressible = (compressibility.value() > SMALL);
3
4rho = thermo.rho();
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());
11surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU));
12volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p_rgh));
13
14surfaceScalarField phig(-rhorAUf*ghf*fvc::snGrad(rho)*mesh.magSf());
15
16surfaceScalarField phiHbyA
17(
18 "phiHbyA",
19 (
20 fvc::flux(rho*HbyA)
21 + MRF.zeroFilter(rhorAUf*fvc::ddtCorr(rho, U, phi))
22 )
23 + phig
24);
25
26MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
27
28// Update the pressure BCs to ensure flux consistency
30
31fvc::makeRelative(phiHbyA, rho, U);
32
33fvScalarMatrix p_rghDDtEqn
34(
35 fvc::ddt(rho) + psi*correction(fvm::ddt(p_rgh))
36 + fvc::div(phiHbyA)
37 ==
38 fvOptions(psi, p_rgh, rho.name())
39);
40
41while (pimple.correctNonOrthogonal())
42{
43 fvScalarMatrix p_rghEqn
44 (
46 - fvm::laplacian(rhorAUf, p_rgh)
47 );
48
49 p_rghEqn.setReference
50 (
52 compressible ? getRefCellValue(p_rgh, pRefCell) : pRefValue
53 );
54
55 p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
56
57 if (pimple.finalNonOrthogonalIter())
58 {
59 // Calculate the conservative fluxes
60 phi = phiHbyA + p_rghEqn.flux();
61
62 // Explicitly relax pressure for momentum corrector
63 p_rgh.relax();
64
65 // Correct the momentum source with the pressure gradient flux
66 // calculated from the relaxed pressure
67 U = HbyA + rAU*fvc::reconstruct((phig + p_rghEqn.flux())/rhorAUf);
68 U.correctBoundaryConditions();
69 fvOptions.correct(U);
70 K = 0.5*magSqr(U);
71 }
72}
73
75
77
78if (p_rgh.needReference())
79{
80 if (!compressible)
81 {
82 p += dimensionedScalar
83 (
84 "p",
85 p.dimensions(),
86 pRefValue - getRefCellValue(p, pRefCell)
87 );
88 }
89 else
90 {
91 p += (initialMass - fvc::domainIntegrate(psi*p))
93 thermo.correctRho(psi*p - psip0, rhoMin, rhoMax);
94 rho = thermo.rho();
95 p_rgh = p - rho*gh;
96 p_rgh.correctBoundaryConditions();
97 }
98}
99else
100{
101 thermo.correctRho(psi*p - psip0, rhoMin, rhoMax);
102}
103
104#include "rhoEqn.H"
105#include "compressibleContinuityErrs.H"
106
107rho = thermo.rho();
108
109// Correct rhoUf if the mesh is moving
110fvc::correctRhoUf(rhoUf, rho, U, phi);
111
112if (thermo.dpdt())
113{
114 dpdt = fvc::ddt(p);
115
116 if (mesh.moving())
117 {
118 dpdt -= fvc::div(fvc::meshPhi(rho, U), p);
119 }
120}
CGAL::Exact_predicates_exact_constructions_kernel K
volScalarField & p_rgh
fv::options & fvOptions
surfaceScalarField & phi
const scalar pRefValue
const surfaceScalarField & ghf
const label pRefCell
IOMRFZoneList & MRF
const volScalarField & gh
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
rhoUf
Definition: pEqn.H:89
volScalarField & p
const volScalarField & psi
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
fvVectorMatrix & UEqn
Definition: UEqn.H:13
phiHbyA
Definition: pcEqn.H:73
HbyA
Definition: pcEqn.H:74
const volScalarField psip0(psi *p)
fvScalarMatrix p_rghDDtEqn(fvc::ddt(rho)+psi *correction(fvm::ddt(p_rgh))+fvc::div(phiHbyA)==fvOptions(psi, p_rgh, rho.name()))
dynamicFvMesh & mesh
const dimensionedScalar rhoMin
const dimensionedScalar rhoMax
volScalarField & dpdt
const pressureControl & pressureControl
dimensionedScalar compressibility
Definition: pEqn.H:1
bool compressible
Definition: pEqn.H:2
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
dimensionedScalar initialMass
Definition: createFields.H:57