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/heatTransfer/buoyantBoussinesqSimpleFoam/pEqn.H
2019-09-26 11:21:23 +01:00

54 lines
1.3 KiB
C++

{
p.boundaryField().updateCoeffs();
volScalarField rAU("rAU", 1.0/UEqn().A());
surfaceScalarField rAUf("(1|A(U))", fvc::interpolate(rAU));
U = rAU*UEqn().H();
UEqn.clear();
simple.calcSteadyConsistentFlux(phi, U);
adjustPhi(phi, U, p);
surfaceScalarField buoyancyPhi =
rAUf*fvc::interpolate(rhok)*(g & mesh.Sf());
phi += buoyancyPhi;
while (simple.correctNonOrthogonal())
{
fvScalarMatrix pEqn
(
fvm::laplacian
(
rAUf/simple.aCoeff(U.name()),
p,
"laplacian(rAU," + p.name() + ')'
)
==
fvc::div(phi)
);
pEqn.setReference(pRefCell, pRefValue);
pEqn.solve();
if (simple.finalNonOrthogonalIter())
{
// Calculate the conservative fluxes
phi -= pEqn.flux();
}
// U += rAU*fvc::reconstruct((buoyancyPhi - pEqn.flux())/rAUf);
// U.correctBoundaryConditions();
}
# include "continuityErrs.H"
// Explicitly relax pressure for momentum corrector
p.relax();
// Correct the momentum source with the pressure gradient flux
// calculated from the relaxed pressure
U += rAU*(fvc::reconstruct(buoyancyPhi/rAUf) - fvc::grad(p));
U.correctBoundaryConditions();
}