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/buoyantSimpleFoam/pEqn.H

58 lines
1.6 KiB
C++
Raw Normal View History

{
2010-08-26 14:22:03 +00:00
rho = thermo.rho();
volScalarField rUA = 1.0/UEqn().A();
surfaceScalarField rhorUAf("(rho*(1|A(U)))", fvc::interpolate(rho*rUA));
U = rUA*UEqn().H();
UEqn.clear();
phi = fvc::interpolate(rho)*(fvc::interpolate(U) & mesh.Sf());
bool closedVolume = adjustPhi(phi, U, p);
surfaceScalarField buoyancyPhi =
rhorUAf*fvc::interpolate(rho)*(g & mesh.Sf());
phi += buoyancyPhi;
while (simple.correctNonOrthogonal())
{
2010-08-26 14:22:03 +00:00
fvScalarMatrix pEqn
(
fvm::laplacian(rhorUAf, p) == fvc::div(phi)
);
2010-08-26 14:22:03 +00:00
pEqn.setReference(pRefCell, pRefValue);
pEqn.solve();
if (simple.finalNonOrthogonalIter())
2010-08-26 14:22:03 +00:00
{
// For closed-volume cases adjust the pressure and density levels
// to obey overall mass continuity
if (closedVolume)
{
p += (initialMass - fvc::domainIntegrate(psi*p))
/fvc::domainIntegrate(psi);
}
2010-08-26 14:22:03 +00:00
// Calculate the conservative fluxes
phi -= pEqn.flux();
2010-08-26 14:22:03 +00:00
// Explicitly relax pressure for momentum corrector
p.relax();
2010-08-26 14:22:03 +00:00
// Correct the momentum source with the pressure gradient flux
// calculated from the relaxed pressure
U += rUA*fvc::reconstruct((buoyancyPhi - pEqn.flux())/rhorUAf);
U.correctBoundaryConditions();
}
}
#include "continuityErrs.H"
rho = thermo.rho();
rho.relax();
Info<< "rho max/min : " << max(rho).value() << " " << min(rho).value()
<< endl;
}