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

62 lines
1.3 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();
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
{
fvScalarMatrix p_rghEqn
(
fvm::laplacian(rAUf, p_rgh) == fvc::div(phi)
);
p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
p_rghEqn.solve
(
mesh.solver
(
p_rgh.select(corr == nCorr-1 && nonOrth == nNonOrthCorr)
)
);
if (nonOrth == nNonOrthCorr)
{
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;
}
}