rhoPimpleFoam: Updated transonic option to be consistent with sonicFoam

Improves stability on start-up and allows running at slightly larger time-steps.
This commit is contained in:
Henry Weller
2017-05-12 11:04:38 +01:00
committed by Andrew Heather
parent 8f1a4f792d
commit 7b825d0817
6 changed files with 145 additions and 121 deletions

View File

@ -1,3 +1,5 @@
rho = thermo.rho();
volScalarField rAU(1.0/UEqn.A()); volScalarField rAU(1.0/UEqn.A());
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU)); surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU));
volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p)); volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
@ -10,7 +12,7 @@ if (pimple.nCorrPISO() <= 1)
surfaceScalarField phiHbyA surfaceScalarField phiHbyA
( (
"phiHbyA", "phiHbyA",
fvc::flux(rho*HbyA) fvc::interpolate(rho)*fvc::flux(HbyA)
+ rhorAUf*fvc::ddtCorr(rho, U, phi) + rhorAUf*fvc::ddtCorr(rho, U, phi)
); );
@ -26,7 +28,8 @@ if (pimple.transonic())
"phid", "phid",
(fvc::interpolate(psi)/fvc::interpolate(rho))*phiHbyA (fvc::interpolate(psi)/fvc::interpolate(rho))*phiHbyA
); );
phiHbyA -= fvc::interpolate(p)*phid;
phiHbyA -= fvc::interpolate(psi*p)*phiHbyA/fvc::interpolate(rho);
while (pimple.correctNonOrthogonal()) while (pimple.correctNonOrthogonal())
{ {
@ -84,9 +87,11 @@ U.correctBoundaryConditions();
fvOptions.correct(U); fvOptions.correct(U);
K = 0.5*magSqr(U); K = 0.5*magSqr(U);
pressureControl.limit(p); if (pressureControl.limit(p))
p.correctBoundaryConditions(); {
rho = thermo.rho(); p.correctBoundaryConditions();
rho = thermo.rho();
}
if (!pimple.transonic()) if (!pimple.transonic())
{ {

View File

@ -1,3 +1,5 @@
rho = thermo.rho();
volScalarField rAU(1.0/UEqn.A()); volScalarField rAU(1.0/UEqn.A());
volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1())); volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()));
volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p)); volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
@ -11,7 +13,7 @@ surfaceScalarField phiHbyA
( (
"phiHbyA", "phiHbyA",
( (
fvc::flux(rho*HbyA) fvc::interpolate(rho)*fvc::flux(HbyA)
+ fvc::interpolate(rho*rAU)*fvc::ddtCorr(rho, U, phi) + fvc::interpolate(rho*rAU)*fvc::ddtCorr(rho, U, phi)
) )
); );
@ -33,7 +35,7 @@ if (pimple.transonic())
phiHbyA += phiHbyA +=
fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf() fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf()
- fvc::interpolate(p)*phid; - fvc::interpolate(psi*p)*phiHbyA/fvc::interpolate(rho);
HbyA -= (rAU - rAtU)*fvc::grad(p); HbyA -= (rAU - rAtU)*fvc::grad(p);
@ -96,9 +98,11 @@ U.correctBoundaryConditions();
fvOptions.correct(U); fvOptions.correct(U);
K = 0.5*magSqr(U); K = 0.5*magSqr(U);
pressureControl.limit(p); if (pressureControl.limit(p))
p.correctBoundaryConditions(); {
rho = thermo.rho(); p.correctBoundaryConditions();
rho = thermo.rho();
}
if (!pimple.transonic()) if (!pimple.transonic())
{ {

View File

@ -1,109 +1,110 @@
volScalarField rAU(1.0/UEqn.A());
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU));
volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
tUEqn.clear();
bool closedVolume = false;
surfaceScalarField phiHbyA("phiHbyA", fvc::interpolate(rho)*fvc::flux(HbyA));
MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
// Update the pressure BCs to ensure flux consistency
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF);
if (simple.transonic())
{ {
volScalarField rAU(1.0/UEqn.A()); surfaceScalarField phid
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU)); (
volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p)); "phid",
tUEqn.clear(); (fvc::interpolate(psi)/fvc::interpolate(rho))*phiHbyA
);
bool closedVolume = false; phiHbyA -= fvc::interpolate(psi*p)*phiHbyA/fvc::interpolate(rho);
surfaceScalarField phiHbyA("phiHbyA", fvc::flux(rho*HbyA)); while (simple.correctNonOrthogonal())
MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
// Update the pressure BCs to ensure flux consistency
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF);
if (simple.transonic())
{ {
surfaceScalarField phid fvScalarMatrix pEqn
( (
"phid", fvc::div(phiHbyA)
(fvc::interpolate(psi)/fvc::interpolate(rho))*phiHbyA + fvm::div(phid, p)
- fvm::laplacian(rhorAUf, p)
==
fvOptions(psi, p, rho.name())
); );
phiHbyA -= fvc::interpolate(p)*phid;
while (simple.correctNonOrthogonal()) // Relax the pressure equation to ensure diagonal-dominance
pEqn.relax();
pEqn.setReference
(
pressureControl.refCell(),
pressureControl.refValue()
);
pEqn.solve();
if (simple.finalNonOrthogonalIter())
{ {
fvScalarMatrix pEqn phi = phiHbyA + pEqn.flux();
(
fvc::div(phiHbyA)
+ fvm::div(phid, p)
- fvm::laplacian(rhorAUf, p)
==
fvOptions(psi, p, rho.name())
);
// Relax the pressure equation to ensure diagonal-dominance
pEqn.relax();
pEqn.setReference
(
pressureControl.refCell(),
pressureControl.refValue()
);
pEqn.solve();
if (simple.finalNonOrthogonalIter())
{
phi = phiHbyA + pEqn.flux();
}
} }
} }
else }
{ else
closedVolume = adjustPhi(phiHbyA, U, p); {
closedVolume = adjustPhi(phiHbyA, U, p);
while (simple.correctNonOrthogonal())
{ while (simple.correctNonOrthogonal())
fvScalarMatrix pEqn {
( fvScalarMatrix pEqn
fvc::div(phiHbyA) (
- fvm::laplacian(rhorAUf, p) fvc::div(phiHbyA)
== - fvm::laplacian(rhorAUf, p)
fvOptions(psi, p, rho.name()) ==
); fvOptions(psi, p, rho.name())
);
pEqn.setReference
( pEqn.setReference
pressureControl.refCell(), (
pressureControl.refValue() pressureControl.refCell(),
); pressureControl.refValue()
);
pEqn.solve();
pEqn.solve();
if (simple.finalNonOrthogonalIter())
{ if (simple.finalNonOrthogonalIter())
phi = phiHbyA + pEqn.flux(); {
} phi = phiHbyA + pEqn.flux();
} }
} }
}
#include "incompressible/continuityErrs.H" #include "incompressible/continuityErrs.H"
// Explicitly relax pressure for momentum corrector // Explicitly relax pressure for momentum corrector
p.relax(); p.relax();
U = HbyA - rAU*fvc::grad(p); U = HbyA - rAU*fvc::grad(p);
U.correctBoundaryConditions(); U.correctBoundaryConditions();
fvOptions.correct(U); fvOptions.correct(U);
pressureControl.limit(p); bool pLimited = pressureControl.limit(p);
// For closed-volume cases adjust the pressure and density levels // For closed-volume cases adjust the pressure and density levels
// to obey overall mass continuity // to obey overall mass continuity
if (closedVolume) if (closedVolume)
{ {
p += (initialMass - fvc::domainIntegrate(psi*p)) p += (initialMass - fvc::domainIntegrate(psi*p))
/fvc::domainIntegrate(psi); /fvc::domainIntegrate(psi);
} }
p.correctBoundaryConditions(); if (pLimited || closedVolume)
{
rho = thermo.rho(); p.correctBoundaryConditions();
}
if (!simple.transonic())
{ rho = thermo.rho();
rho.relax();
} if (!simple.transonic())
{
rho.relax();
} }

View File

@ -1,3 +1,5 @@
rho = thermo.rho();
volScalarField rAU(1.0/UEqn.A()); volScalarField rAU(1.0/UEqn.A());
volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1())); volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()));
volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p)); volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
@ -5,7 +7,7 @@ tUEqn.clear();
bool closedVolume = false; bool closedVolume = false;
surfaceScalarField phiHbyA("phiHbyA", fvc::flux(rho*HbyA)); surfaceScalarField phiHbyA("phiHbyA", fvc::interpolate(rho)*fvc::flux(HbyA));
MRF.makeRelative(fvc::interpolate(rho), phiHbyA); MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
volScalarField rhorAtU("rhorAtU", rho*rAtU); volScalarField rhorAtU("rhorAtU", rho*rAtU);
@ -23,7 +25,7 @@ if (simple.transonic())
phiHbyA += phiHbyA +=
fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf() fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf()
- fvc::interpolate(p)*phid; - fvc::interpolate(psi*p)*phiHbyA/fvc::interpolate(rho);
HbyA -= (rAU - rAtU)*fvc::grad(p); HbyA -= (rAU - rAtU)*fvc::grad(p);
@ -98,7 +100,7 @@ U = HbyA - rAtU*fvc::grad(p);
U.correctBoundaryConditions(); U.correctBoundaryConditions();
fvOptions.correct(U); fvOptions.correct(U);
pressureControl.limit(p); bool pLimited = pressureControl.limit(p);
// For closed-volume cases adjust the pressure and density levels // For closed-volume cases adjust the pressure and density levels
// to obey overall mass continuity // to obey overall mass continuity
@ -108,9 +110,11 @@ if (closedVolume)
/fvc::domainIntegrate(psi); /fvc::domainIntegrate(psi);
} }
p.correctBoundaryConditions(); if (pLimited || closedVolume)
{
p.correctBoundaryConditions();
}
// Recalculate density from the relaxed pressure
rho = thermo.rho(); rho = thermo.rho();
if (!simple.transonic()) if (!simple.transonic())

View File

@ -207,14 +207,24 @@ Foam::pressureControl::pressureControl
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
void Foam::pressureControl::limit(volScalarField& p) const bool Foam::pressureControl::limit(volScalarField& p) const
{ {
Info<< "pressureControl: p max/min " scalar pMax = max(p).value();
<< max(p).value() << " " scalar pMin = min(p).value();
<< min(p).value() << endl;
p = max(p, pMin_); if (pMax > pMax_.value() || pMin < pMin_.value())
p = min(p, pMax_); {
Info<< "pressureControl: p max/min " << pMax << " " << pMin << endl;
p = max(p, pMin_);
p = min(p, pMax_);
return true;
}
else
{
return false;
}
} }

View File

@ -89,8 +89,8 @@ public:
//- Return the pressure reference level //- Return the pressure reference level
inline scalar refValue() const; inline scalar refValue() const;
//- Limit the pressure //- Limit the pressure if necessary and return true if so
void limit(volScalarField& p) const; bool limit(volScalarField& p) const;
}; };