pEqn.H
Go to the documentation of this file.
1{
2 volScalarField rAU("rAU", 1.0/UEqn.A());
3 surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU));
4 volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p_rgh));
5 surfaceScalarField phiHbyA
6 (
7 "phiHbyA",
8 fvc::flux(HbyA)
9 + MRF.zeroFilter(fvc::interpolate(rho*rAU)*fvc::ddtCorr(U, Uf))
10 );
11 MRF.makeRelative(phiHbyA);
12
13 surfaceScalarField phig
14 (
15 (
16 mixture.surfaceTensionForce()
17 - ghf*fvc::snGrad(rho)
18 )*rAUf*mesh.magSf()
19 );
20
22
23 // Update the pressure BCs to ensure flux consistency
25
26 // Make the fluxes relative to the mesh motion
27 fvc::makeRelative(phiHbyA, U);
28
29 tmp<fvScalarMatrix> p_rghEqnComp1;
30 tmp<fvScalarMatrix> p_rghEqnComp2;
31
32 if (pimple.transonic())
33 {
34 #include "rhofs.H"
35
36 surfaceScalarField phid1("phid1", fvc::interpolate(psi1)*phi);
37 surfaceScalarField phid2("phid2", fvc::interpolate(psi2)*phi);
38
40 (
41 (
42 fvc::ddt(alpha1, rho1) + fvc::div(alphaPhi1*rho1f)
43 - (fvOptions(alpha1, mixture.thermo1().rho())&rho1)
44 )/rho1
45 - fvc::ddt(alpha1) - fvc::div(alphaPhi1)
46 + (alpha1/rho1)
47 *correction
48 (
49 psi1*fvm::ddt(p_rgh)
50 + fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
51 )
52 );
53 p_rghEqnComp1.ref().relax();
54
56 (
57 (
58 fvc::ddt(alpha2, rho2) + fvc::div(alphaPhi2*rho2f)
59 - (fvOptions(alpha2, mixture.thermo2().rho())&rho2)
60 )/rho2
61 - fvc::ddt(alpha2) - fvc::div(alphaPhi2)
62 + (alpha2/rho2)
63 *correction
64 (
65 psi2*fvm::ddt(p_rgh)
66 + fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
67 )
68 );
69 p_rghEqnComp2.ref().relax();
70 }
71 else
72 {
73 #include "rhofs.H"
74
76 pos(alpha1)
77 *(
78 (
79 fvc::ddt(alpha1, rho1) + fvc::div(alphaPhi1*rho1f)
80 - (fvOptions(alpha1, mixture.thermo1().rho())&rho1)
81 )/rho1
82 - fvc::ddt(alpha1) - fvc::div(alphaPhi1)
83 + (alpha1*psi1/rho1)*correction(fvm::ddt(p_rgh))
84 );
85
87 pos(alpha2)
88 *(
89 (
90 fvc::ddt(alpha2, rho2) + fvc::div(alphaPhi2*rho2f)
91 - (fvOptions(alpha2, mixture.thermo2().rho())&rho2)
92 )/rho2
93 - fvc::ddt(alpha2) - fvc::div(alphaPhi2)
94 + (alpha2*psi2/rho2)*correction(fvm::ddt(p_rgh))
95 );
96 }
97
98 if (mesh.moving())
99 {
100 p_rghEqnComp1.ref() += fvc::div(mesh.phi())*alpha1;
101 p_rghEqnComp2.ref() += fvc::div(mesh.phi())*alpha2;
102 }
103
104 p_rghEqnComp1.ref() *= pos(alpha1);
105 p_rghEqnComp2.ref() *= pos(alpha2);
106
107 if (pimple.transonic())
108 {
109 p_rghEqnComp1.ref().relax();
110 p_rghEqnComp2.ref().relax();
111 }
112
113 // Cache p_rgh prior to solve for density update
114 volScalarField p_rgh_0(p_rgh);
115
116 while (pimple.correctNonOrthogonal())
117 {
118 fvScalarMatrix p_rghEqnIncomp
119 (
120 fvc::div(phiHbyA)
121 - fvm::laplacian(rAUf, p_rgh)
122 );
123
124 solve
125 (
126 p_rghEqnComp1() + p_rghEqnComp2() + p_rghEqnIncomp,
127 mesh.solver(p_rgh.select(pimple.finalInnerIter()))
128 );
129
130 if (pimple.finalNonOrthogonalIter())
131 {
132 p = max(p_rgh + (alpha1*rho1 + alpha2*rho2)*gh, pMin);
133 p_rgh = p - (alpha1*rho1 + alpha2*rho2)*gh;
134
135 dgdt =
136 (
139 );
140
141 phi = phiHbyA + p_rghEqnIncomp.flux();
142
143 U = HbyA
144 + rAU*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
145 U.correctBoundaryConditions();
146 fvOptions.correct(U);
147 }
148 }
149
150 // Correct Uf if the mesh is moving
151 {
152 Uf = fvc::interpolate(U);
153 surfaceVectorField n(mesh.Sf()/mesh.magSf());
154 Uf += n*(fvc::absolute(phi, U)/mesh.magSf() - (n & Uf));
155 }
156
157
158 // Update densities from change in p_rgh
159 mixture.thermo1().correctRho(psi1*(p_rgh - p_rgh_0));
160 mixture.thermo2().correctRho(psi2*(p_rgh - p_rgh_0));
161
163
164 // Correct p_rgh for consistency with p and the updated densities
166 p_rgh.correctBoundaryConditions();
167
168 K = 0.5*magSqr(U);
169}
CGAL::Exact_predicates_exact_constructions_kernel K
label n
Y[inertIndex] max(0.0)
volScalarField & p_rgh
fv::options & fvOptions
surfaceScalarField & phi
const surfaceScalarField & ghf
IOMRFZoneList & MRF
const volScalarField & gh
pimpleControl & pimple
const volScalarField & alpha1
const volScalarField & psi2
volScalarField & rho2
const volScalarField & alpha2
const surfaceScalarField & alphaPhi1
const volScalarField & psi1
const dimensionedScalar & pMin
volScalarField & rho1
const surfaceScalarField & alphaPhi2
U
Definition: pEqn.H:72
volScalarField & p
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
fvVectorMatrix & UEqn
Definition: UEqn.H:13
phiHbyA
Definition: pcEqn.H:73
HbyA
Definition: pcEqn.H:74
surfaceScalarField rho2f(fvc::interpolate(rho2))
surfaceScalarField rho1f(fvc::interpolate(rho1))
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
tmp< fvScalarMatrix > p_rghEqnComp1
Definition: pEqn.H:29
volScalarField p_rgh_0(p_rgh)
tmp< fvScalarMatrix > p_rghEqnComp2
Definition: pEqn.H:30
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
CEqn solve()
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39