Loading...
Searching...
No Matches
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());
11mesh.interpolate(rAU);
12
13surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU));
14volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p_rgh));
15
16surfaceScalarField phig(-rhorAUf*ghf*fvc::snGrad(rho)*mesh.magSf());
17
18surfaceScalarField phiHbyA
19(
20 "phiHbyA",
21 fvc::flux(rho*HbyA) + phig
22);
23
25{
26 fvc::makeRelative(phiHbyA,rho, U);
27 oversetAdjustPhi(phiHbyA, U);
28 fvc::makeAbsolute(phiHbyA,rho, U);
29}
30
31
32MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
33
34// Update the pressure BCs to ensure flux consistency
36
37fvScalarMatrix p_rghDDtEqn
38(
39 fvc::ddt(rho) + psi*correction(fvm::ddt(p_rgh))
40 + fvc::div(phiHbyA)
41 ==
42 fvOptions(psi, p_rgh, rho.name())
43);
44
45while (pimple.correctNonOrthogonal())
46{
47 fvScalarMatrix p_rghEqn
48 (
50 - fvm::laplacian(rhorAUf, p_rgh)
51 );
52
53 p_rghEqn.solve(p_rgh.select(pimple.finalInnerIter()));
54
55 if (pimple.finalNonOrthogonalIter())
56 {
57 // Calculate the conservative fluxes
58 phi = phiHbyA + p_rghEqn.flux();
59
60 // Explicitly relax pressure for momentum corrector
61 p_rgh.relax();
62
63 // Correct the momentum source with the pressure gradient flux
64 // calculated from the relaxed pressure
65 U =
66 cellMask*
67 (
68 HbyA + rAU*fvc::reconstruct((phig + p_rghEqn.flux())/rhorAUf)
69 );
70 U.correctBoundaryConditions();
71 fvOptions.correct(U);
72 K = 0.5*magSqr(U);
73 }
74}
75
77
78#include "rhoEqn.H"
79#include "compressibleContinuityErrs.H"
80
81if (p_rgh.needReference())
82{
83 if (!compressible)
84 {
85 p += dimensionedScalar
86 (
87 "p",
88 p.dimensions(),
89 pRefValue - getRefCellValue(p, pRefCell)
90 );
91 }
92 else
93 {
94 p += (initialMass - fvc::domainIntegrate(psi*p))
96 thermo.correctRho(psi*p - psip0);
97 rho = thermo.rho();
98 p_rgh = p - rho*gh;
99 }
100}
101else
102{
103 thermo.correctRho(psi*p - psip0);
104}
105
106rho = thermo.rho();
107
108{
109 fvc::correctRhoUf(rhoUf, rho, U, phi);
110}
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}
121
122phi *= faceMask;
CGAL::Exact_predicates_exact_constructions_kernel K
volScalarField & p_rgh
fv::options & fvOptions
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
volScalarField & p
const volScalarField & psi
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
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
autoPtr< surfaceVectorField > rhoUf
volScalarField & dpdt
dimensionedScalar compressibility
Definition pEqn.H:1
bool compressible
Definition pEqn.H:2
tmp< volScalarField > rAU
faceMask
Definition setCellMask.H:43
dimensionedScalar initialMass