Files
openfoam/applications/solvers/multiphase/compressibleTwoPhaseEulerFoam/pEqn.H

158 lines
4.0 KiB
C

{
rho1 = rho10 + psi1*p;
rho2 = rho20 + psi2*p;
mrfZones.absoluteFlux(phi1.oldTime());
mrfZones.absoluteFlux(phi1);
mrfZones.absoluteFlux(phi2.oldTime());
mrfZones.absoluteFlux(phi2);
surfaceScalarField alpha1f(fvc::interpolate(alpha1));
surfaceScalarField alpha2f(scalar(1) - alpha1f);
volScalarField rAU1(1.0/U1Eqn.A());
volScalarField rAU2(1.0/U2Eqn.A());
surfaceScalarField rAlphaAU1f = fvc::interpolate(alpha1*rAU1);
surfaceScalarField rAlphaAU2f = fvc::interpolate(alpha2*rAU2);
volVectorField HbyA1("HbyA1", U1);
HbyA1 = rAU1*U1Eqn.H();
volVectorField HbyA2("HbyA2", U2);
HbyA2 = rAU2*U2Eqn.H();
surfaceScalarField phiHbyA1
(
"phiHbyA1",
(fvc::interpolate(HbyA1) & mesh.Sf())
+ fvc::ddtPhiCorr(rAU1, alpha1, U1, phi1)
+ fvc::interpolate((1.0/rho1)*rAU1*dragCoeff)*phi2
+ rAlphaAU1f*(g & mesh.Sf())
);
mrfZones.relativeFlux(phiHbyA1);
surfaceScalarField phiHbyA2
(
"phiHbyA2",
(fvc::interpolate(HbyA2) & mesh.Sf())
+ fvc::ddtPhiCorr(rAU2, alpha2, U2, phi2)
+ fvc::interpolate((1.0/rho2)*rAU2*dragCoeff)*phi1
+ rAlphaAU2f*(g & mesh.Sf())
);
mrfZones.relativeFlux(phiHbyA2);
mrfZones.relativeFlux(phi1.oldTime());
mrfZones.relativeFlux(phi1);
mrfZones.relativeFlux(phi2.oldTime());
mrfZones.relativeFlux(phi2);
surfaceScalarField phiHbyA("phiHbyA", alpha1f*phiHbyA1 + alpha2f*phiHbyA2);
HbyA1 += (1.0/rho1)*rAU1*dragCoeff*U2;
HbyA2 += (1.0/rho2)*rAU2*dragCoeff*U1;
surfaceScalarField Dp
(
"Dp",
mag
(
alpha1f*rAlphaAU1f/fvc::interpolate(rho1)
+ alpha2f*rAlphaAU2f/fvc::interpolate(rho2)
)
);
tmp<fvScalarMatrix> pEqnComp1;
tmp<fvScalarMatrix> pEqnComp2;
//if (transonic)
//{
//}
//else
{
surfaceScalarField phid1("phid1", fvc::interpolate(psi1)*phi1);
surfaceScalarField phid2("phid2", fvc::interpolate(psi2)*phi2);
pEqnComp1 =
fvc::ddt(rho1) + psi1*correction(fvm::ddt(p))
+ fvc::div(phid1, p)
- fvc::Sp(fvc::div(phid1), p);
pEqnComp2 =
fvc::ddt(rho2) + psi2*correction(fvm::ddt(p))
+ fvc::div(phid2, p)
- fvc::Sp(fvc::div(phid2), p);
}
while (pimple.correctNonOrthogonal())
{
fvScalarMatrix pEqnIncomp
(
fvc::div(phiHbyA)
- fvm::laplacian(Dp, p)
);
solve
(
(
(alpha1/rho1)*pEqnComp1()
+ (alpha2/rho2)*pEqnComp2()
)
+ pEqnIncomp,
mesh.solver(p.select(pimple.finalInnerIter()))
);
if (pimple.finalNonOrthogonalIter())
{
surfaceScalarField mSfGradp = pEqnIncomp.flux()/Dp;
phi1 = phiHbyA1 + rAlphaAU1f*mSfGradp/fvc::interpolate(rho1);
phi2 = phiHbyA2 + rAlphaAU2f*mSfGradp/fvc::interpolate(rho2);
phi = alpha1f*phi1 + alpha2f*phi2;
dgdt =
(
pos(alpha2)*(pEqnComp2 & p)/rho2
- pos(alpha1)*(pEqnComp1 & p)/rho1
);
p.relax();
mSfGradp = pEqnIncomp.flux()/Dp;
U1 = HbyA1
+ fvc::reconstruct
(
rAlphaAU1f
*(
(g & mesh.Sf())
+ mSfGradp/fvc::interpolate(rho1)
)
);
U1.correctBoundaryConditions();
U2 = HbyA2
+ fvc::reconstruct
(
rAlphaAU2f
*(
(g & mesh.Sf())
+ mSfGradp/fvc::interpolate(rho2)
)
);
U2.correctBoundaryConditions();
U = alpha1*U1 + alpha2*U2;
}
}
p = max(p, pMin);
rho1 = rho10 + psi1*p;
rho2 = rho20 + psi2*p;
K1 = 0.5*magSqr(U1);
K2 = 0.5*magSqr(U1);
dpdt = fvc::ddt(p);
}