mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
simpleFoam/SRFSimpleFoam: Added support for SIMPLEC
SIMPLEC (SIMPLE-consistent) is selected by setting "consistent" option true/yes:
SIMPLE
{
nNonOrthogonalCorrectors 0;
consistent yes;
}
which relaxes the pressure in a "consistent" manner and additional
relaxation of the pressure is not generally necessary. In addition
convergence of the p-U system is better and reliable with less
aggressive relaxation of the momentum equation, e.g. for the motorbike
tutorial:
relaxationFactors
{
equations
{
U 0.9;
k 0.7;
omega 0.7;
}
}
The cost per iteration is marginally higher but the convergence rate is
better so the number of iterations can be reduced.
The SIMPLEC algorithm also provides benefit for cases with large
body-forces, e.g. SRF, see tutorials/incompressible/SRFSimpleFoam/mixer
and feature request http://www.openfoam.org/mantisbt/view.php?id=1714
This commit is contained in:
@ -2,17 +2,28 @@
|
||||
volScalarField rAUrel(1.0/UrelEqn().A());
|
||||
volVectorField HbyA("HbyA", Urel);
|
||||
HbyA = rAUrel*UrelEqn().H();
|
||||
UrelEqn.clear();
|
||||
|
||||
surfaceScalarField phiHbyA("phiHbyA", fvc::interpolate(HbyA) & mesh.Sf());
|
||||
adjustPhi(phiHbyA, Urel, p);
|
||||
|
||||
tmp<volScalarField> rAtUrel(rAUrel);
|
||||
|
||||
if (simple.consistent())
|
||||
{
|
||||
rAtUrel = 1.0/(1.0/rAUrel - UrelEqn().H1());
|
||||
phiHbyA +=
|
||||
fvc::interpolate(rAtUrel() - rAUrel)*fvc::snGrad(p)*mesh.magSf();
|
||||
HbyA -= (rAUrel - rAtUrel())*fvc::grad(p);
|
||||
}
|
||||
|
||||
UrelEqn.clear();
|
||||
|
||||
// Non-orthogonal pressure corrector loop
|
||||
while (simple.correctNonOrthogonal())
|
||||
{
|
||||
fvScalarMatrix pEqn
|
||||
(
|
||||
fvm::laplacian(rAUrel, p) == fvc::div(phiHbyA)
|
||||
fvm::laplacian(rAtUrel(), p) == fvc::div(phiHbyA)
|
||||
);
|
||||
|
||||
pEqn.setReference(pRefCell, pRefValue);
|
||||
@ -31,7 +42,7 @@
|
||||
p.relax();
|
||||
|
||||
// Momentum corrector
|
||||
Urel = HbyA - rAUrel*fvc::grad(p);
|
||||
Urel = HbyA - rAtUrel()*fvc::grad(p);
|
||||
Urel.correctBoundaryConditions();
|
||||
fvOptions.correct(Urel);
|
||||
}
|
||||
|
||||
@ -2,20 +2,29 @@
|
||||
volScalarField rAU(1.0/UEqn().A());
|
||||
volVectorField HbyA("HbyA", U);
|
||||
HbyA = rAU*UEqn().H();
|
||||
UEqn.clear();
|
||||
|
||||
surfaceScalarField phiHbyA("phiHbyA", fvc::interpolate(HbyA) & mesh.Sf());
|
||||
|
||||
fvOptions.makeRelative(phiHbyA);
|
||||
|
||||
adjustPhi(phiHbyA, U, p);
|
||||
|
||||
tmp<volScalarField> rAtU(rAU);
|
||||
|
||||
if (simple.consistent())
|
||||
{
|
||||
rAtU = 1.0/(1.0/rAU - UEqn().H1());
|
||||
phiHbyA +=
|
||||
fvc::interpolate(rAtU() - rAU)*fvc::snGrad(p)*mesh.magSf();
|
||||
HbyA -= (rAU - rAtU())*fvc::grad(p);
|
||||
}
|
||||
|
||||
UEqn.clear();
|
||||
|
||||
// Non-orthogonal pressure corrector loop
|
||||
while (simple.correctNonOrthogonal())
|
||||
{
|
||||
fvScalarMatrix pEqn
|
||||
(
|
||||
fvm::laplacian(rAU, p) == fvc::div(phiHbyA)
|
||||
fvm::laplacian(rAtU(), p) == fvc::div(phiHbyA)
|
||||
);
|
||||
|
||||
pEqn.setReference(pRefCell, pRefValue);
|
||||
@ -34,7 +43,7 @@
|
||||
p.relax();
|
||||
|
||||
// Momentum corrector
|
||||
U = HbyA - rAU*fvc::grad(p);
|
||||
U = HbyA - rAtU()*fvc::grad(p);
|
||||
U.correctBoundaryConditions();
|
||||
fvOptions.correct(U);
|
||||
}
|
||||
|
||||
Reference in New Issue
Block a user