/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | foam-extend: Open Source CFD \\ / O peration | Version: 4.0 \\ / A nd | Web: http://www.foam-extend.org \\/ M anipulation | For copyright notice see file Copyright ------------------------------------------------------------------------------- License This file is part of foam-extend. 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 Free Software Foundation, either version 3 of the License, or (at your option) any later version. 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 along with foam-extend. If not, see . Application MRFSimpleFoam Description Steady-state solver for incompressible, turbulent flow of non-Newtonian fluids with MRF regions. \*---------------------------------------------------------------------------*/ #include "fvCFD.H" #include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H" #include "incompressible/RAS/RASModel/RASModel.H" #include "MRFZones.H" #include "simpleControl.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[]) { # include "setRootCase.H" # include "createTime.H" # include "createMesh.H" simpleControl simple(mesh); # include "createFields.H" # include "initContinuityErrs.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // Info<< "\nStarting time loop\n" << endl; while (simple.loop()) { Info<< "Time = " << runTime.timeName() << nl << endl; // Pressure-velocity SIMPLE corrector { // Momentum predictor tmp UEqn ( fvm::div(phi, U) + turbulence->divDevReff() ); mrfZones.addCoriolis(UEqn()); UEqn().relax(); solve(UEqn() == -fvc::grad(p)); volScalarField rAU = 1.0/UEqn().A(); U = rAU*UEqn().H(); UEqn.clear(); phi = fvc::interpolate(U, "interpolate(HbyA)") & mesh.Sf(); mrfZones.relativeFlux(phi); adjustPhi(phi, U, p); // Non-orthogonal pressure corrector loop while (simple.correctNonOrthogonal()) { fvScalarMatrix pEqn ( fvm::laplacian(rAU, p) == fvc::div(phi) ); pEqn.setReference(pRefCell, pRefValue); pEqn.solve(); if (simple.finalNonOrthogonalIter()) { phi -= pEqn.flux(); } } # include "continuityErrs.H" // Explicitly relax pressure for momentum corrector p.relax(); // Momentum corrector U -= rAU*fvc::grad(p); U.correctBoundaryConditions(); } turbulence->correct(); runTime.write(); Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s" << " ClockTime = " << runTime.elapsedClockTime() << " s" << nl << endl; } Info<< "End\n" << endl; return(0); } // ************************************************************************* //