/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2017 OpenFOAM Foundation Copyright (C) 2019 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 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] \ | As required by user selection \endplaintable Note The motion frequency of this solver can be influenced by the presence of "updateControl" and "updateInterval" in the dynamicMeshDict. \*---------------------------------------------------------------------------*/ #include "gpufvCFD.H" #include "dynamicFvMesh.H" #include "singlePhaseTransportModel.H" #include "turbulentTransportModel.H" #include "pimpleControl.H" #include "gpuCorrectPhi.H" #include "gpufvOptions.H" #include "gpulocalEulerDdtScheme.H" #include "gpufvcSmooth.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[]) { argList::addNote ( "Transient solver for incompressible, turbulent flow" " of Newtonian fluids on a moving mesh." ); #include "gpupostProcess.H" #include "addCheckCaseOptions.H" #include "setRootCaseLists.H" #include "createTime.H" #include "createDynamicFvMesh.H" Foam::gpufvMesh devicemesh(mesh); #include "initContinuityErrs.H" #include "createDyMControls.H" #include "createFields.H" #include "gpucreateUfIfPresent.H" #include "gpuCourantNo.H" #include "setInitialDeltaT.H" turbulence->validate(); if (!LTS) { #include "gpuCourantNo.H" #include "setInitialDeltaT.H" } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // Info<< "\nStarting time loop\n" << endl; while (runTime.run()) { #include "readDyMControls.H" if (LTS) { #include "setRDeltaT.H" } else { #include "gpuCourantNo.H" #include "setDeltaT.H" } ++runTime; Info<< "Time = " << runTime.timeName() << nl << endl; // --- Pressure-velocity PIMPLE corrector loop while (pimple.loop()) { // todo in loop if (pimple.firstIter() || moveMeshOuterCorrectors) { // Do any mesh changes mesh.controlledUpdate(); if (mesh.changing()) { //update mesh geometry devicemesh.updateGpuGeom(); U.correctBoundaryConditions(); MRF.update(); if (correctPhi) { // Calculate absolute flux // from the mapped surface velocity phi = devicemesh.Sf() & Uf(); #include "correctPhi.H" // Make the flux relative to the mesh motion fvc::makeRelative(phi, U); } if (checkMeshCourantNo) { #include "gpumeshCourantNo.H" //todo } } } #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; } // ************************************************************************* //