Files
openfoam/applications/solvers/compressible/rhoSimpleFoam/pcEqn.H
Henry Weller 77b03e2e0c Specialized dotInterpolate for the efficient calculation of flux fields
e.g. (fvc::interpolate(HbyA) & mesh.Sf()) -> fvc::flux(HbyA)

This removes the need to create an intermediate face-vector field when
computing fluxes which is more efficient, reduces the peak storage and
improved cache coherency in addition to providing a simpler and cleaner
API.
2016-04-06 20:20:53 +01:00

119 lines
2.7 KiB
C

volScalarField rAU(1.0/UEqn.A());
volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()));
volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
tUEqn.clear();
bool closedVolume = false;
if (simple.transonic())
{
surfaceScalarField phid
(
"phid",
fvc::interpolate(psi)
*fvc::flux(HbyA)
);
MRF.makeRelative(fvc::interpolate(psi), phid);
surfaceScalarField phic
(
"phic",
fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf()
);
HbyA -= (rAU - rAtU)*fvc::grad(p);
volScalarField rhorAtU("rhorAtU", rho*rAtU);
while (simple.correctNonOrthogonal())
{
fvScalarMatrix pEqn
(
fvm::div(phid, p)
+ fvc::div(phic)
- fvm::laplacian(rhorAtU, p)
==
fvOptions(psi, p, rho.name())
);
// Relax the pressure equation to maintain diagonal dominance
pEqn.relax();
pEqn.setReference(pRefCell, pRefValue);
pEqn.solve();
if (simple.finalNonOrthogonalIter())
{
phi == phic + pEqn.flux();
}
}
}
else
{
surfaceScalarField phiHbyA("phiHbyA", fvc::flux(rho*HbyA));
MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
closedVolume = adjustPhi(phiHbyA, U, p);
phiHbyA += fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf();
HbyA -= (rAU - rAtU)*fvc::grad(p);
volScalarField rhorAtU("rhorAtU", rho*rAtU);
// Update the pressure BCs to ensure flux consistency
constrainPressure(p, rho, U, phiHbyA, rhorAtU, MRF);
while (simple.correctNonOrthogonal())
{
fvScalarMatrix pEqn
(
fvc::div(phiHbyA)
- fvm::laplacian(rhorAtU, p)
==
fvOptions(psi, p, rho.name())
);
pEqn.setReference(pRefCell, pRefValue);
pEqn.solve();
if (simple.finalNonOrthogonalIter())
{
phi = phiHbyA + pEqn.flux();
}
}
}
// The incompressibe form of the continuity error check is appropriate for
// steady-state compressible also.
#include "incompressible/continuityErrs.H"
// Explicitly relax pressure for momentum corrector
p.relax();
U = HbyA - rAtU*fvc::grad(p);
U.correctBoundaryConditions();
fvOptions.correct(U);
// 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);
}
// Recalculate density from the relaxed pressure
rho = thermo.rho();
rho = max(rho, rhoMin);
rho = min(rho, rhoMax);
if (!simple.transonic())
{
rho.relax();
}
Info<< "rho max/min : " << max(rho).value() << " " << min(rho).value() << endl;