- relocate some standard functionality to TimePaths to allow a lighter means of managing time directories without using the entire Time mechanism. - optional enableLibs for Time construction (default is on) and a corresponding argList::noLibs() and "-no-libs" option STYLE: - mark Time::outputTime() as deprecated MAY-2016 - use pre-increment for runTime, although there is no difference in behaviour or performance.
168 lines
4.6 KiB
C
168 lines
4.6 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
|
|
\\/ 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 3 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, see <http://www.gnu.org/licenses/>.
|
|
|
|
Application
|
|
pimpleFoam.C
|
|
|
|
Group
|
|
grpIncompressibleSolvers
|
|
|
|
Description
|
|
Transient solver for incompressible, turbulent flow of Newtonian fluids
|
|
on a moving mesh.
|
|
|
|
\heading Solver details
|
|
The solver uses the PIMPLE (merged PISO-SIMPLE) algorithm to solve the
|
|
continuity equation:
|
|
|
|
\f[
|
|
\div \vec{U} = 0
|
|
\f]
|
|
|
|
and momentum equation:
|
|
|
|
\f[
|
|
\ddt{\vec{U}} + \div \left( \vec{U} \vec{U} \right) - \div \gvec{R}
|
|
= - \grad p + \vec{S}_U
|
|
\f]
|
|
|
|
Where:
|
|
\vartable
|
|
\vec{U} | Velocity
|
|
p | Pressure
|
|
\vec{R} | Stress tensor
|
|
\vec{S}_U | Momentum source
|
|
\endvartable
|
|
|
|
Sub-models include:
|
|
- turbulence modelling, i.e. laminar, RAS or LES
|
|
- run-time selectable MRF and finite volume options, e.g. explicit porosity
|
|
|
|
\heading Required fields
|
|
\plaintable
|
|
U | Velocity [m/s]
|
|
p | Kinematic pressure, p/rho [m2/s2]
|
|
\<turbulence fields\> | As required by user selection
|
|
\endplaintable
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "fvCFD.H"
|
|
#include "dynamicFvMesh.H"
|
|
#include "singlePhaseTransportModel.H"
|
|
#include "turbulentTransportModel.H"
|
|
#include "pimpleControl.H"
|
|
#include "CorrectPhi.H"
|
|
#include "fvOptions.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
int main(int argc, char *argv[])
|
|
{
|
|
#include "postProcess.H"
|
|
|
|
#include "addCheckCaseOptions.H"
|
|
#include "setRootCase.H"
|
|
#include "createTime.H"
|
|
#include "createDynamicFvMesh.H"
|
|
#include "initContinuityErrs.H"
|
|
#include "createDyMControls.H"
|
|
#include "createFields.H"
|
|
#include "createUfIfPresent.H"
|
|
#include "CourantNo.H"
|
|
#include "setInitialDeltaT.H"
|
|
|
|
turbulence->validate();
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
Info<< "\nStarting time loop\n" << endl;
|
|
|
|
while (runTime.run())
|
|
{
|
|
#include "readDyMControls.H"
|
|
#include "CourantNo.H"
|
|
#include "setDeltaT.H"
|
|
|
|
++runTime;
|
|
|
|
Info<< "Time = " << runTime.timeName() << nl << endl;
|
|
|
|
// --- Pressure-velocity PIMPLE corrector loop
|
|
while (pimple.loop())
|
|
{
|
|
if (pimple.firstIter() || moveMeshOuterCorrectors)
|
|
{
|
|
mesh.update();
|
|
|
|
if (mesh.changing())
|
|
{
|
|
MRF.update();
|
|
|
|
if (correctPhi)
|
|
{
|
|
// Calculate absolute flux
|
|
// from the mapped surface velocity
|
|
phi = mesh.Sf() & Uf();
|
|
|
|
#include "correctPhi.H"
|
|
|
|
// Make the flux relative to the mesh motion
|
|
fvc::makeRelative(phi, U);
|
|
}
|
|
|
|
if (checkMeshCourantNo)
|
|
{
|
|
#include "meshCourantNo.H"
|
|
}
|
|
}
|
|
}
|
|
|
|
#include "UEqn.H"
|
|
|
|
// --- Pressure corrector loop
|
|
while (pimple.correct())
|
|
{
|
|
#include "pEqn.H"
|
|
}
|
|
|
|
if (pimple.turbCorr())
|
|
{
|
|
laminarTransport.correct();
|
|
turbulence->correct();
|
|
}
|
|
}
|
|
|
|
runTime.write();
|
|
|
|
runTime.printExecutionTime(Info);
|
|
}
|
|
|
|
Info<< "End\n" << endl;
|
|
|
|
return 0;
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|