-
Notifications
You must be signed in to change notification settings - Fork 5
/
pEqn.H
143 lines (123 loc) · 4.24 KB
/
pEqn.H
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
{
Psi1 = rho_gas/p_rgh;
Psi2 = rho_foam/p_rgh;
volScalarField rAU("rAU", 1.0/UEqn.A());
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU));
volVectorField HbyA("HbyA", U);
HbyA = rAU*UEqn.H();
surfaceScalarField phiHbyA
(
"phiHbyA",
(fvc::interpolate(HbyA) & mesh.Sf())
+ fvc::interpolate(rho*rAU)*fvc::ddtCorr(U, phi)
);
surfaceScalarField phig
(
(
fvc::interpolate(interface.sigmaK())*fvc::snGrad(alpha1)
- ghf*fvc::snGrad(rho)
)*rAUf*mesh.magSf()
);
phiHbyA += phig;
// Update the fixedFluxPressure BCs to ensure flux consistency
setSnGrad<fixedFluxPressureFvPatchScalarField>
(
p_rgh.boundaryField(),
(
phiHbyA.boundaryField()
- (mesh.Sf().boundaryField() & U.boundaryField())
)/(mesh.magSf().boundaryField()*rAUf.boundaryField())
);
tmp<fvScalarMatrix> p_rghEqnComp1;
tmp<fvScalarMatrix> p_rghEqnComp2;
if (pimple.transonic())
{
surfaceScalarField phid1("phid1", fvc::interpolate(Psi1)*phi);
surfaceScalarField phid2("phid2", fvc::interpolate(Psi2)*phi);
p_rghEqnComp1 =
fvc::ddt(rho_gas) + fvc::div(phi, rho_gas)
- fvc::Sp(fvc::div(phi), rho_gas)
+ correction
(
Psi1*fvm::ddt(p_rgh)
+ fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
);
deleteDemandDrivenData(p_rghEqnComp1().faceFluxCorrectionPtr());
p_rghEqnComp1().relax();
p_rghEqnComp2 =
fvc::ddt(rho_foam) + fvc::div(phi, rho_foam)
- fvc::Sp(fvc::div(phi), rho_foam)
+ correction
(
Psi2*fvm::ddt(p_rgh)
+ fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
);
deleteDemandDrivenData(p_rghEqnComp2().faceFluxCorrectionPtr());
p_rghEqnComp2().relax();
}
else
{
p_rghEqnComp1 =
(
fvc::ddt(alpha1, rho_gas)
+ fvc::div(fvc::interpolate(alpha1)*phi, rho_gas)
- fvc::Sp(fvc::ddt(alpha1)
+ fvc::div(fvc::interpolate(alpha1)*phi), rho_gas)
)/rho_gas
+ (alpha1*Psi1/rho_gas)*correction(fvm::ddt(p_rgh));
p_rghEqnComp2 =
(
fvc::ddt(alpha2, rho_foam)
+ fvc::div(fvc::interpolate(alpha2)*phi, rho_foam)
- fvc::Sp(fvc::ddt(alpha2)
+ fvc::div(fvc::interpolate(alpha2)*phi), rho_foam)
)/rho_foam
+ (alpha2*Psi2/rho_foam)*correction(fvm::ddt(p_rgh));
}
// Cache p_rgh prior to solve for density update
volScalarField p_rgh_0(p_rgh);
while (pimple.correctNonOrthogonal())
{
fvScalarMatrix p_rghEqnIncomp
(
fvc::div(phiHbyA)
- fvm::laplacian(rAUf, p_rgh)
);
solve
(
(
// (max(alpha1, scalar(0))/rho_gas)*
p_rghEqnComp1()
+
// (max((1-alpha1), scalar(0))/rho_foam)*
p_rghEqnComp2()
)
+ p_rghEqnIncomp,
mesh.solver(p_rgh.select(pimple.finalInnerIter()))
);
if (pimple.finalNonOrthogonalIter())
{
p = max(p_rgh + (alpha1*rho_gas + (1.0-alpha1)*rho_foam)*gh, pMin);
p_rgh = p - (alpha1*rho_gas + (1.0-alpha1)*rho_foam)*gh;
dimensionedScalar dumdgdt ("dumdgdt",dimensionSet(1,-3,0,0,0,0,0), 1.0);
dgdt =
dumdgdt*
(
pos((1.0-alpha1))*(p_rghEqnComp2 & p_rgh)/rho_foam
- pos(alpha1)*(p_rghEqnComp1 & p_rgh)/rho_gas
);
phi = phiHbyA + p_rghEqnIncomp.flux();
U = HbyA
+ rAU*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
U.correctBoundaryConditions();
}
}
p = max(p_rgh + (alpha1*rho_gas + (1.0-alpha1)*rho_foam)*gh, pMin);
// Update densities from change in p_rgh
rho_gas += Psi1*(p_rgh - p_rgh_0);
rho_foam += Psi2*(p_rgh - p_rgh_0);
rho = alpha1*rho_gas + (1-alpha1)*rho_foam;
K = 0.5*magSqr(U);
Info<< "max(U) " << max(mag(U)).value() << endl;
Info<< "min(p_rgh) " << min(p_rgh).value() << endl;
}