/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2014 OpenFOAM Foundation Copyright (C) 2016-2017 OpenCFD Ltd. ------------------------------------------------------------------------------- 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 . Application overSimpleFoam Group grpIncompressibleSolvers Description Steady-state solver for incompressible flows with turbulence modelling. \heading Solver details The solver uses the SIMPLE algorithm to solve the continuity equation: \f[ \div \vec{U} = 0 \f] and momentum equation: \f[ \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 \heading Required fields \plaintable U | Velocity [m/s] p | Kinematic pressure, p/rho [m2/s2] \ | As required by user selection \endplaintable \*---------------------------------------------------------------------------*/ #include "fvCFD.H" #include "singlePhaseTransportModel.H" #include "turbulentTransportModel.H" #include "simpleControl.H" #include "fvOptions.H" #include "dynamicFvMesh.H" #include "cellCellStencilObject.H" #include "localMin.H" #include "interpolationCellPoint.H" #include "fvMeshSubset.H" #include "transform.H" #include "oversetAdjustPhi.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[]) { argList::addNote ( "Steady-state solver for incompressible, turbulent flow" ); #define CREATE_MESH createUpdatedDynamicFvMesh.H #include "postProcess.H" #include "setRootCaseLists.H" #include "createTime.H" #include "createUpdatedDynamicFvMesh.H" #include "createControl.H" #include "createFields.H" #include "createOversetFields.H" #include "createFvOptions.H" #include "initContinuityErrs.H" turbulence->validate(); // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // Info<< "\nStarting time loop\n" << endl; while (simple.loop()) { Info<< "Time = " << runTime.timeName() << nl << endl; // --- Pressure-velocity SIMPLE corrector { #include "UEqn.H" #include "pEqn.H" } laminarTransport.correct(); turbulence->correct(); runTime.write(); runTime.printExecutionTime(Info); } Info<< "End\n" << endl; return 0; } // ************************************************************************* //