Backported pimpleFoam (solutionControl and setFluxRequired) from OF 3.0.1
This commit is contained in:
parent
99fdf7833b
commit
a704cb308e
4 changed files with 25 additions and 46 deletions
|
@ -7,34 +7,20 @@ tmp<fvVectorMatrix> UEqn
|
|||
+ turbulence->divDevReff(U)
|
||||
);
|
||||
|
||||
if (oCorr == nOuterCorr - 1)
|
||||
{
|
||||
if (mesh.solutionDict().relax("UFinal"))
|
||||
{
|
||||
UEqn().relax(mesh.solutionDict().relaxationFactor("UFinal"));
|
||||
}
|
||||
else
|
||||
{
|
||||
UEqn().relax(1);
|
||||
}
|
||||
}
|
||||
else
|
||||
{
|
||||
UEqn().relax();
|
||||
}
|
||||
UEqn().relax
|
||||
(
|
||||
mesh.solutionDict().relaxationFactor(U.select(pimple.finalIter()))
|
||||
);
|
||||
|
||||
volScalarField rUA = 1.0/UEqn().A();
|
||||
|
||||
if (momentumPredictor)
|
||||
if (pimple.momentumPredictor())
|
||||
{
|
||||
if (oCorr == nOuterCorr-1)
|
||||
{
|
||||
solve(UEqn() == -fvc::grad(p), mesh.solutionDict().solver("UFinal"));
|
||||
}
|
||||
else
|
||||
{
|
||||
solve(UEqn() == -fvc::grad(p));
|
||||
}
|
||||
solve
|
||||
(
|
||||
UEqn() == -fvc::grad(p),
|
||||
mesh.solutionDict().solver((U.select(pimple.finalIter())))
|
||||
);
|
||||
}
|
||||
else
|
||||
{
|
||||
|
|
|
@ -32,7 +32,7 @@ volVectorField U
|
|||
label pRefCell = 0;
|
||||
scalar pRefValue = 0.0;
|
||||
setRefCell(p, mesh.solutionDict().subDict("PIMPLE"), pRefCell, pRefValue);
|
||||
|
||||
mesh.schemesDict().setFluxRequired(p.name());
|
||||
|
||||
singlePhaseTransportModel laminarTransport(U, phi);
|
||||
|
||||
|
|
|
@ -1,6 +1,6 @@
|
|||
U = rUA*UEqn().H();
|
||||
|
||||
if (nCorr <= 1)
|
||||
if (pimple.nCorrPISO() <= 1)
|
||||
{
|
||||
UEqn.clear();
|
||||
}
|
||||
|
@ -11,7 +11,7 @@ phi = (fvc::interpolate(U) & mesh.Sf())
|
|||
adjustPhi(phi, U, p);
|
||||
|
||||
// Non-orthogonal pressure corrector loop
|
||||
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
|
||||
while (pimple.correctNonOrthogonal())
|
||||
{
|
||||
// Pressure corrector
|
||||
fvScalarMatrix pEqn
|
||||
|
@ -20,22 +20,12 @@ for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
|
|||
);
|
||||
|
||||
pEqn.setReference(pRefCell, pRefValue);
|
||||
|
||||
if
|
||||
pEqn.solve
|
||||
(
|
||||
oCorr == nOuterCorr - 1
|
||||
&& corr == nCorr - 1
|
||||
&& nonOrth == nNonOrthCorr
|
||||
)
|
||||
{
|
||||
pEqn.solve(mesh.solutionDict().solver("pFinal"));
|
||||
}
|
||||
else
|
||||
{
|
||||
pEqn.solve();
|
||||
}
|
||||
mesh.solutionDict().solver(p.select(pimple.finalInnerIter()))
|
||||
);
|
||||
|
||||
if (nonOrth == nNonOrthCorr)
|
||||
if (pimple.finalNonOrthogonalIter())
|
||||
{
|
||||
phi -= pEqn.flux();
|
||||
}
|
||||
|
@ -44,7 +34,7 @@ for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
|
|||
#include "continuityErrs.H"
|
||||
|
||||
// Explicitly relax pressure for momentum corrector except for last corrector
|
||||
if (oCorr != nOuterCorr-1)
|
||||
if (!pimple.finalIter())
|
||||
{
|
||||
p.relax();
|
||||
}
|
||||
|
|
|
@ -35,6 +35,7 @@ Description
|
|||
#include "fvCFD.H"
|
||||
#include "singlePhaseTransportModel.H"
|
||||
#include "turbulenceModel.H"
|
||||
#include "pimpleControl.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
|
@ -43,6 +44,9 @@ int main(int argc, char *argv[])
|
|||
# include "setRootCase.H"
|
||||
# include "createTime.H"
|
||||
# include "createMesh.H"
|
||||
|
||||
pimpleControl pimple(mesh);
|
||||
|
||||
# include "createFields.H"
|
||||
# include "initContinuityErrs.H"
|
||||
|
||||
|
@ -51,7 +55,6 @@ int main(int argc, char *argv[])
|
|||
while (runTime.run())
|
||||
{
|
||||
# include "readTimeControls.H"
|
||||
# include "readPIMPLEControls.H"
|
||||
# include "CourantNo.H"
|
||||
# include "setDeltaT.H"
|
||||
|
||||
|
@ -60,9 +63,9 @@ int main(int argc, char *argv[])
|
|||
Info<< "Time = " << runTime.timeName() << nl << endl;
|
||||
|
||||
// --- Pressure-velocity PIMPLE corrector loop
|
||||
for (int oCorr = 0; oCorr < nOuterCorr; oCorr++)
|
||||
while (pimple.loop())
|
||||
{
|
||||
if (nOuterCorr != 1)
|
||||
if (!pimple.firstIter())
|
||||
{
|
||||
p.storePrevIter();
|
||||
}
|
||||
|
@ -70,7 +73,7 @@ int main(int argc, char *argv[])
|
|||
# include "UEqn.H"
|
||||
|
||||
// --- PISO loop
|
||||
for (int corr = 0; corr < nCorr; corr++)
|
||||
while (pimple.correct())
|
||||
{
|
||||
# include "pEqn.H"
|
||||
}
|
||||
|
|
Reference in a new issue