This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/incompressible/pisoFoam/pisoFoam.C

157 lines
4.8 KiB
C++
Raw Normal View History

/*---------------------------------------------------------------------------*\
========= |
2013-12-11 16:09:41 +00:00
\\ / F ield | foam-extend: Open Source CFD
2016-06-20 15:00:40 +00:00
\\ / O peration | Version: 4.0
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
2013-12-11 16:09:41 +00:00
This file is part of foam-extend.
2013-12-11 16:09:41 +00:00
foam-extend is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
2013-12-11 16:09:41 +00:00
Free Software Foundation, either version 3 of the License, or (at your
option) any later version.
2013-12-11 16:09:41 +00:00
foam-extend 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
2013-12-11 16:09:41 +00:00
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
Application
2010-08-26 14:22:03 +00:00
pisoFoam
Description
2017-01-05 08:46:27 +00:00
Transient solver for incompressible, turbulent flow.
2010-08-26 14:22:03 +00:00
Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
2017-01-05 08:46:27 +00:00
Consistent formulation without time-step and relaxation dependence by Jasak
and Tukovic.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
2010-08-26 14:22:03 +00:00
#include "singlePhaseTransportModel.H"
#include "turbulenceModel.H"
#include "pisoControl.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
# include "setRootCase.H"
# include "createTime.H"
# include "createMesh.H"
pisoControl piso(mesh);
# include "createFields.H"
# include "initContinuityErrs.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
2010-08-26 14:22:03 +00:00
while (runTime.loop())
{
Info<< "Time = " << runTime.timeName() << nl << endl;
# include "CourantNo.H"
// Pressure-velocity PISO corrector
{
// Momentum predictor
2017-01-05 08:46:27 +00:00
// Time-derivative matrix
fvVectorMatrix ddtUEqn(fvm::ddt(U));
// Convection-diffusion matrix
fvVectorMatrix HUEqn
(
2017-01-05 08:46:27 +00:00
fvm::div(phi, U)
+ turbulence->divDevReff()
);
if (piso.momentumPredictor())
{
2017-01-05 08:46:27 +00:00
solve(relax(ddtUEqn + HUEqn) == -fvc::grad(p));
}
// --- PISO loop
while (piso.correct())
{
2017-01-05 08:46:27 +00:00
// Prepare clean 1/a_p without time derivative and
// under-relaxation contribution
volScalarField rAU = 1.0/HUEqn.A();
2017-01-05 08:46:27 +00:00
// Calculate U from convection-diffusion matrix
U = rAU*HUEqn.H();
2017-01-05 08:46:27 +00:00
// Consistently calculate flux
piso.calcTransientConsistentFlux(phi, U, rAU, ddtUEqn);
// Global flux balance
adjustPhi(phi, U, p);
// Non-orthogonal pressure corrector loop
while (piso.correctNonOrthogonal())
{
// Pressure corrector
fvScalarMatrix pEqn
(
2017-01-05 08:46:27 +00:00
fvm::laplacian
(
fvc::interpolate(rAU)/piso.aCoeff(U.name()),
2017-01-05 08:46:27 +00:00
p,
"laplacian(rAU," + p.name() + ')'
)
==
fvc::div(phi)
);
pEqn.setReference(pRefCell, pRefValue);
pEqn.solve
(
mesh.solutionDict().solver
(
p.select(piso.finalInnerIter())
)
);
if (piso.finalNonOrthogonalIter())
{
phi -= pEqn.flux();
}
}
# include "continuityErrs.H"
2017-01-05 08:46:27 +00:00
// Consistently reconstruct velocity after pressure
// equation. Note: flux is made relative inside the function
piso.reconstructTransientVelocity(U, phi, ddtUEqn, rAU, p);
}
}
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;
}
// ************************************************************************* //