This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/multiphase/twoLiquidMixingFoam/p_rghEqn.H
2011-10-12 10:32:38 +01:00

56 lines
1.2 KiB
C

{
volScalarField rUA = 1.0/UEqn.A();
surfaceScalarField rUAf = fvc::interpolate(rUA);
U = rUA*UEqn.H();
surfaceScalarField phiU
(
"phiU",
(fvc::interpolate(U) & mesh.Sf())
+ fvc::ddtPhiCorr(rUA, rho, U, phi)
);
phi = phiU - ghf*fvc::snGrad(rho)*rUAf*mesh.magSf();
for(int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
{
fvScalarMatrix p_rghEqn
(
fvm::laplacian(rUAf, p_rgh) == fvc::div(phi)
);
p_rghEqn.setReference(p_rghRefCell, p_rghRefValue);
p_rghEqn.solve
(
mesh.solutionDict().solver
(
p_rgh.name()
+ ((corr == nCorr - 1 && nonOrth == nNonOrthCorr) ? "Final" : "")
)
);
if (nonOrth == nNonOrthCorr)
{
phi -= p_rghEqn.flux();
}
}
p == p_rgh + rho*gh;
if (p_rgh.needReference())
{
p += dimensionedScalar
(
"p",
p.dimensions(),
pRefValue - getRefCellValue(p, p_rghRefCell)
);
}
#include "continuityErrs.H"
U += rUA*fvc::reconstruct((phi - phiU)/rUAf);
U.correctBoundaryConditions();
}