2010-05-12 13:27:55 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
|
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
|
|
\\ / O peration |
|
|
|
|
\\ / A nd | Copyright held by original author
|
|
|
|
\\/ M anipulation |
|
|
|
|
-------------------------------------------------------------------------------
|
|
|
|
License
|
|
|
|
This file is part of OpenFOAM.
|
|
|
|
|
|
|
|
OpenFOAM is free software; you can redistribute it and/or modify it
|
|
|
|
under the terms of the GNU General Public License as published by the
|
|
|
|
Free Software Foundation; either version 2 of the License, or (at your
|
|
|
|
option) any later version.
|
|
|
|
|
|
|
|
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
|
|
|
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
|
|
|
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
|
|
|
for more details.
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OpenFOAM; if not, write to the Free Software Foundation,
|
|
|
|
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
|
|
|
|
|
|
|
Application
|
2010-08-26 14:22:03 +00:00
|
|
|
pisoFoam
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
Description
|
2010-08-26 14:22:03 +00:00
|
|
|
Transient solver for incompressible flow.
|
|
|
|
|
|
|
|
Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#include "fvCFD.H"
|
2010-08-26 14:22:03 +00:00
|
|
|
#include "singlePhaseTransportModel.H"
|
|
|
|
#include "turbulenceModel.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
int main(int argc, char *argv[])
|
|
|
|
{
|
|
|
|
# include "setRootCase.H"
|
|
|
|
|
|
|
|
# include "createTime.H"
|
|
|
|
# include "createMesh.H"
|
|
|
|
# include "createFields.H"
|
|
|
|
# include "initContinuityErrs.H"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
Info<< "\nStarting time loop\n" << endl;
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
while (runTime.loop())
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
Info<< "Time = " << runTime.timeName() << nl << endl;
|
|
|
|
|
|
|
|
# include "readPISOControls.H"
|
|
|
|
# include "CourantNo.H"
|
|
|
|
|
|
|
|
// Pressure-velocity PISO corrector
|
|
|
|
{
|
|
|
|
// Momentum predictor
|
|
|
|
|
|
|
|
fvVectorMatrix UEqn
|
|
|
|
(
|
|
|
|
fvm::ddt(U)
|
|
|
|
+ fvm::div(phi, U)
|
|
|
|
+ turbulence->divDevReff(U)
|
|
|
|
);
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
UEqn.relax();
|
|
|
|
|
2010-05-12 13:27:55 +00:00
|
|
|
if (momentumPredictor)
|
|
|
|
{
|
|
|
|
solve(UEqn == -fvc::grad(p));
|
|
|
|
}
|
|
|
|
|
|
|
|
// --- PISO loop
|
|
|
|
|
|
|
|
for (int corr=0; corr<nCorr; corr++)
|
|
|
|
{
|
|
|
|
volScalarField rUA = 1.0/UEqn.A();
|
|
|
|
|
|
|
|
U = rUA*UEqn.H();
|
2010-08-26 14:22:03 +00:00
|
|
|
phi = (fvc::interpolate(U) & mesh.Sf())
|
2010-05-12 13:27:55 +00:00
|
|
|
+ fvc::ddtPhiCorr(rUA, U, phi);
|
|
|
|
|
|
|
|
adjustPhi(phi, U, p);
|
|
|
|
|
|
|
|
// Non-orthogonal pressure corrector loop
|
|
|
|
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
|
|
|
|
{
|
|
|
|
// Pressure corrector
|
|
|
|
|
|
|
|
fvScalarMatrix pEqn
|
|
|
|
(
|
|
|
|
fvm::laplacian(rUA, p) == fvc::div(phi)
|
|
|
|
);
|
|
|
|
|
|
|
|
pEqn.setReference(pRefCell, pRefValue);
|
2010-08-26 14:22:03 +00:00
|
|
|
|
|
|
|
if
|
|
|
|
(
|
|
|
|
corr == nCorr-1
|
|
|
|
&& nonOrth == nNonOrthCorr
|
|
|
|
)
|
|
|
|
{
|
2011-08-14 16:39:59 +00:00
|
|
|
pEqn.solve(mesh.solutionDict().solver("pFinal"));
|
2010-08-26 14:22:03 +00:00
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
pEqn.solve();
|
|
|
|
}
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
if (nonOrth == nNonOrthCorr)
|
|
|
|
{
|
|
|
|
phi -= pEqn.flux();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
# include "continuityErrs.H"
|
|
|
|
|
|
|
|
U -= rUA*fvc::grad(p);
|
|
|
|
U.correctBoundaryConditions();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
turbulence->correct();
|
|
|
|
|
|
|
|
runTime.write();
|
|
|
|
|
|
|
|
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
|
|
|
|
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
|
|
|
|
<< nl << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< "End\n" << endl;
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
return 0;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|