Loading...
Searching...
No Matches
pEqn.H
Go to the documentation of this file.
1rho = thermo.rho();
2rho.clamp_range(rhoMin, rhoMax);
3rho.relax();
4
5volScalarField rAU(1.0/UEqn.A());
6surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU));
7volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
8
9if (pimple.nCorrPISO() <= 1)
10{
11 tUEqn.clear();
12}
13
14if (pimple.transonic())
15{
16 surfaceScalarField phid
17 (
18 "phid",
19 fvc::interpolate(psi)
20 *(
21 fvc::flux(HbyA)
22 + MRF.zeroFilter
23 (
24 rhorAUf*fvc::ddtCorr(rho, U, rhoUf)/fvc::interpolate(rho)
25 )
26 )
27 );
28
29 fvc::makeRelative(phid, psi, U);
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(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 fvc::flux(rho*HbyA)
58 + MRF.zeroFilter(rhorAUf*fvc::ddtCorr(rho, U, rhoUf))
59 );
60
61 fvc::makeRelative(phiHbyA, rho, U);
62 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
63
64 // Update the pressure BCs to ensure flux consistency
66
67 while (pimple.correctNonOrthogonal())
68 {
69 fvScalarMatrix pEqn
70 (
71 fvm::ddt(psi, p)
72 + fvc::div(phiHbyA)
73 - fvm::laplacian(rhorAUf, p)
74 ==
75 parcels.Srho()
76 + fvOptions(psi, p, rho.name())
77 );
78
79 pEqn.solve(p.select(pimple.finalInnerIter()));
80
81 if (pimple.finalNonOrthogonalIter())
82 {
83 phi = phiHbyA + pEqn.flux();
84 }
85 }
86}
87
88#include "rhoEqn.H"
89#include "compressibleContinuityErrs.H"
90
91// Explicitly relax pressure for momentum corrector
92p.relax();
93
94// Recalculate density from the relaxed pressure
95rho = thermo.rho();
96rho.clamp_range(rhoMin, rhoMax);
97rho.relax();
98Info<< "rho min/max : " << min(rho).value() << " " << max(rho).value() << endl;
99
100U = HbyA - rAU*fvc::grad(p);
101U.correctBoundaryConditions();
102fvOptions.correct(U);
103K = 0.5*magSqr(U);
104
105{
106 rhoUf = fvc::interpolate(rho*U);
107 surfaceVectorField n(mesh.Sf()/mesh.magSf());
108 rhoUf += n*(fvc::absolute(phi, rho, U)/mesh.magSf() - (n & rhoUf));
109}
110
111if (thermo.dpdt())
112{
113 dpdt = fvc::ddt(p);
114
115 if (mesh.moving())
116 {
117 dpdt -= fvc::div(fvc::meshPhi(rho, U), p);
118 }
119}
CGAL::Exact_predicates_exact_constructions_kernel K
label n
fv::options & fvOptions
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
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
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
autoPtr< surfaceVectorField > rhoUf
const dimensionedScalar rhoMin
const dimensionedScalar rhoMax
volScalarField & dpdt
tmp< volScalarField > rAU
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)