ENH: Updated bubbleFoam solver and tutorial

This commit is contained in:
andy 2011-04-18 15:05:30 +01:00
parent 6ab6476be7
commit 39b30e2429
5 changed files with 47 additions and 29 deletions

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2004-2010 OpenCFD Ltd.
\\ / A nd | Copyright (C) 2004-2011 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -34,6 +34,7 @@ Description
#include "nearWallDist.H"
#include "wallFvPatch.H"
#include "Switch.H"
#include "pimpleControl.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -47,6 +48,8 @@ int main(int argc, char *argv[])
#include "createFields.H"
#include "initContinuityErrs.H"
pimpleControl pimple(mesh);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
@ -58,26 +61,32 @@ int main(int argc, char *argv[])
#include "readBubbleFoamControls.H"
#include "CourantNo.H"
#include "alphaEqn.H"
#include "liftDragCoeffs.H"
#include "UEqns.H"
// --- PISO loop
for (int corr=0; corr<nCorr; corr++)
// --- Pressure-velocity PIMPLE corrector loop
for (pimple.start(); pimple.loop(); pimple++)
{
#include "pEqn.H"
#include "alphaEqn.H"
#include "liftDragCoeffs.H"
#include "UEqns.H"
if (correctAlpha)
// --- PISO loop
for (int corr=0; corr<pimple.nCorr(); corr++)
{
#include "alphaEqn.H"
#include "pEqn.H"
if (correctAlpha)
{
#include "alphaEqn.H"
}
}
#include "DDtU.H"
if (pimple.turbCorr())
{
#include "kEpsilon.H"
}
}
#include "DDtU.H"
#include "kEpsilon.H"
#include "write.H"
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"

View File

@ -192,4 +192,4 @@
label pRefCell = 0;
scalar pRefValue = 0.0;
setRefCell(p, mesh.solutionDict().subDict("PISO"), pRefCell, pRefValue);
setRefCell(p, mesh.solutionDict().subDict("PIMPLE"), pRefCell, pRefValue);

View File

@ -42,7 +42,7 @@
alphaf*rUaAf/rhoa + betaf*rUbAf/rhob
);
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
for (int nonOrth=0; nonOrth<=pimple.nNonOrthCorr(); nonOrth++)
{
fvScalarMatrix pEqn
(
@ -50,9 +50,13 @@
);
pEqn.setReference(pRefCell, pRefValue);
pEqn.solve();
if (nonOrth == nNonOrthCorr)
pEqn.solve
(
mesh.solver(p.select(pimple.finalInnerIter(corr, nonOrth)))
);
if (nonOrth == pimple.nNonOrthCorr())
{
surfaceScalarField SfGradp(pEqn.flux()/Dp);

View File

@ -1,5 +1,4 @@
#include "readPISOControls.H"
int nAlphaCorr(readInt(pisoDict.lookup("nAlphaCorr")));
int nAlphaCorr(readInt(pimple.dict().lookup("nAlphaCorr")));
Switch correctAlpha(pisoDict.lookup("correctAlpha"));
Switch correctAlpha(pimple.dict().lookup("correctAlpha"));

View File

@ -22,27 +22,33 @@ solvers
solver PCG;
preconditioner DIC;
tolerance 1e-10;
relTol 0.1;
}
pFinal
{
$p;
tolerance 1e-10;
relTol 0;
}
"(Ua|Ub|k|epsilon)"
alpha
{
solver PBiCG;
preconditioner DILU;
tolerance 1e-05;
relTol 0;
tolerance 1e-10;
relTol 0.1;
}
"(alpha|beta)"
alphaFinal
{
solver PBiCG;
preconditioner DILU;
$alpha;
tolerance 1e-10;
relTol 0;
}
}
PISO
PIMPLE
{
nCorrectors 2;
nNonOrthogonalCorrectors 0;