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

70 lines
1.4 KiB
C++
Raw Permalink Normal View History

2014-12-16 22:37:18 +00:00
{
U = UEqn.H()/UEqn.A();
# include "limitU.H"
while (pimple.correctNonOrthogonal())
{
// Calculate phi for boundary conditions
phi = rhof*
2014-12-16 22:37:18 +00:00
(
(fvc::interpolate(U) & mesh.Sf())
- fvc::meshPhi(rho, U)
);
surfaceScalarField phid2 = rhoReff/rhof*phi;
surfaceScalarField phid("phid", psisf/rhof*phi);
2014-12-16 22:37:18 +00:00
// Store pressure for under-relaxation
p.storePrevIter();
volScalarField divPhid
(
"divPhid",
fvc::div(phid)
);
2014-12-16 22:37:18 +00:00
fvScalarMatrix pEqn
(
fvm::ddt(psis, p)
2014-12-16 22:37:18 +00:00
+ fvm::div(phid, p)
// Convective flux relaxation terms
+ fvm::SuSp(-divPhid, p)
+ divPhid*p
+ fvc::div(phid2)
- fvm::laplacian(rho*rUA, p)
2014-12-16 22:37:18 +00:00
);
pEqn.solve
2014-12-16 22:37:18 +00:00
(
mesh.solutionDict().solver(p.select(pimple.finalInnerIter()))
);
2014-12-16 22:37:18 +00:00
// Calculate the flux
if (pimple.finalNonOrthogonalIter())
2014-12-16 22:37:18 +00:00
{
phi = phid2 + pEqn.flux();
2014-12-16 22:37:18 +00:00
}
// Bound the pressure
if (min(p) < pMin || max(p) > pMax)
{
p.max(pMin);
p.min(pMax);
p.correctBoundaryConditions();
}
// Relax the pressure
p.relax();
}
# include "compressibleContinuityErrs.H"
U -= fvc::grad(p)/UEqn.A();
U.correctBoundaryConditions();
# include "limitU.H"
}