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

63 lines
1.4 KiB
C

{
volScalarField rAU(1.0/UEqn.A());
surfaceScalarField rAUf(fvc::interpolate(rAU));
U = rAU*UEqn.H();
surfaceScalarField phiU
(
"phiU",
(fvc::interpolate(U) & mesh.Sf())
+ fvc::ddtPhiCorr(rAU, rho, U, phi)
);
adjustPhi(phiU, U, p_rgh);
phi =
phiU
+ (
fvc::interpolate(interface.sigmaK())*fvc::snGrad(alpha1)
- ghf*fvc::snGrad(rho)
)*rAUf*mesh.magSf();
Pair<tmp<volScalarField> > vDotP = twoPhaseProperties->vDotP();
const volScalarField& vDotcP = vDotP[0]();
const volScalarField& vDotvP = vDotP[1]();
while (pimple.correctNonOrthogonal())
{
fvScalarMatrix p_rghEqn
(
fvc::div(phi) - fvm::laplacian(rAUf, p_rgh)
- (vDotvP - vDotcP)*(pSat - rho*gh) + fvm::Sp(vDotvP - vDotcP, p_rgh)
);
p_rghEqn.setReference(pRefCell, pRefValue);
p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
if (pimple.finalNonOrthogonalIter())
{
phi += p_rghEqn.flux();
}
}
U += rAU*fvc::reconstruct((phi - phiU)/rAUf);
U.correctBoundaryConditions();
#include "continuityErrs.H"
p == p_rgh + rho*gh;
if (p_rgh.needReference())
{
p += dimensionedScalar
(
"p",
p.dimensions(),
pRefValue - getRefCellValue(p, pRefCell)
);
p_rgh = p - rho*gh;
}
}