openfoam/applications/solvers/multiphase/interPhaseChangeFoam/overInterPhaseChangeDyMFoam/pEqn.H
Mark Olesen ef44df91f2 ENH: support direct lookup of solver controls
OLD:
        pEqn.solve(mesh.solver(p.select(piso.finalInnerIter())));
        pEqn.solve(mesh.solver("Yi"));

    NEW:
        pEqn.solve(p.select(piso.finalInnerIter()));
        pEqn.solve("Yi");
2023-12-07 17:42:24 +01:00

92 lines
2.2 KiB
C

{
rAU = 1.0/UEqn.A();
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU));
volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p_rgh));
surfaceScalarField phiHbyA
(
"phiHbyA",
fvc::flux(HbyA)
);
if (p_rgh.needReference())
{
fvc::makeRelative(phiHbyA, U);
adjustPhi(phiHbyA, U, p_rgh);
fvc::makeAbsolute(phiHbyA, U);
}
surfaceScalarField phig
(
(
interface.surfaceTensionForce()
- ghf*fvc::snGrad(rho)
)*faceMask*rAUf*mesh.magSf()
);
phiHbyA += phig;
// Update the pressure BCs to ensure flux consistency
constrainPressure(p_rgh, U, phiHbyA, rAUf);
Pair<tmp<volScalarField>> vDotP = mixture->vDotP();
const volScalarField& vDotcP = vDotP[0]();
const volScalarField& vDotvP = vDotP[1]();
while (pimple.correctNonOrthogonal())
{
fvScalarMatrix p_rghEqn
(
fvc::div(phiHbyA) - fvm::laplacian(rAUf, p_rgh)
- (vDotvP - vDotcP)*(mixture->pSat() - rho*gh)
+ fvm::Sp(vDotvP - vDotcP, p_rgh)
);
//p_rghEqn.setReference(pRefCell, pRefValue);
p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
p_rghEqn.solve(p_rgh.select(pimple.finalInnerIter()));
if (pimple.finalNonOrthogonalIter())
{
phi = phiHbyA + p_rghEqn.flux();
p_rgh.relax();
U =
cellMask
*(HbyA + rAU*fvc::reconstruct((phig + p_rghEqn.flux())/rAUf));
U.correctBoundaryConditions();
fvOptions.correct(U);
}
}
#include "continuityErrs.H"
{
Uf = fvc::interpolate(U);
surfaceVectorField n(mesh.Sf()/mesh.magSf());
Uf += n*(phi/mesh.magSf() - (n & Uf));
}
// Make the fluxes relative to the mesh motion
fvc::makeRelative(phi, U);
// Zero faces H-I for transport Eq after pEq
phi *= faceMask;
p == p_rgh + rho*gh;
if (p_rgh.needReference())
{
p += dimensionedScalar
(
"p",
p.dimensions(),
pRefValue - getRefCellValue(p, pRefCell)
);
p_rgh = p - rho*gh;
}
}