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

104 lines
2.3 KiB
C
Raw Normal View History

2010-08-26 14:22:03 +00:00
rho = thermo.rho();
rho = max(rho, rhoMin);
rho = min(rho, rhoMax);
rho.relax();
volScalarField rUA = 1.0/UEqn().A();
U = rUA*UEqn().H();
UEqn.clear();
2010-08-26 14:22:03 +00:00
bool closedVolume = false;
if (transonic)
{
2010-08-26 14:22:03 +00:00
surfaceScalarField phid
(
2010-08-26 14:22:03 +00:00
"phid",
fvc::interpolate(psi)*(fvc::interpolate(U) & mesh.Sf())
);
2010-08-26 14:22:03 +00:00
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
{
2010-08-26 14:22:03 +00:00
fvScalarMatrix pEqn
(
fvm::div(phid, p)
- fvm::laplacian(rho*rUA, p)
);
// Relax the pressure equation to ensure diagonal-dominance
pEqn.relax(mesh.solutionDict().relaxationFactor("pEqn"));
2010-08-26 14:22:03 +00:00
pEqn.setReference(pRefCell, pRefValue);
// retain the residual from the first iteration
if (nonOrth == 0)
{
eqnResidual = pEqn.solve().initialResidual();
maxResidual = max(eqnResidual, maxResidual);
}
else
{
pEqn.solve();
}
if (nonOrth == nNonOrthCorr)
{
phi == pEqn.flux();
}
}
2010-08-26 14:22:03 +00:00
}
else
{
phi = fvc::interpolate(rho)*(fvc::interpolate(U) & mesh.Sf());
closedVolume = adjustPhi(phi, U, p);
2010-08-26 14:22:03 +00:00
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
{
2010-08-26 14:22:03 +00:00
fvScalarMatrix pEqn
(
fvm::laplacian(rho*rUA, p) == fvc::div(phi)
);
pEqn.setReference(pRefCell, pRefValue);
// Retain the residual from the first iteration
if (nonOrth == 0)
{
eqnResidual = pEqn.solve().initialResidual();
maxResidual = max(eqnResidual, maxResidual);
}
else
{
pEqn.solve();
}
if (nonOrth == nNonOrthCorr)
{
phi -= pEqn.flux();
}
}
}
2010-08-26 14:22:03 +00:00
#include "incompressible/continuityErrs.H"
// Explicitly relax pressure for momentum corrector
p.relax();
2010-08-26 14:22:03 +00:00
U -= rUA*fvc::grad(p);
U.correctBoundaryConditions();
// For closed-volume cases adjust the pressure and density levels
// to obey overall mass continuity
if (closedVolume)
{
2010-08-26 14:22:03 +00:00
p += (initialMass - fvc::domainIntegrate(psi*p))
/fvc::domainIntegrate(psi);
}
2010-08-26 14:22:03 +00:00
rho = thermo.rho();
rho = max(rho, rhoMin);
rho = min(rho, rhoMax);
rho.relax();
Info<< "rho max/min : " << max(rho).value() << " " << min(rho).value() << endl;