initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / tutorials / incompressible / MRFSimpleFoam / MRFSimpleFoam / MRFSimpleFoam.C
blob9225fe4cb7c20e4456307054f6ade7afa7aec3d3
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     MRFSimpleFoam
28 Description
29     Steady-state solver for incompressible, turbulent flow of non-Newtonian
30     fluids with MRF regions.
32 \*---------------------------------------------------------------------------*/
34 #include "fvCFD.H"
35 #include "singlePhaseTransportModel.H"
36 #include "RASModel.H"
37 #include "MRFZones.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 int main(int argc, char *argv[])
43     #include "setRootCase.H"
45     #include "createTime.H"
46     #include "createMesh.H"
47     #include "createFields.H"
48     #include "initContinuityErrs.H"
50     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52     Info<< "\nStarting time loop\n" << endl;
54     while (runTime.loop())
55     {
56         Info<< "Time = " << runTime.timeName() << nl << endl;
58         #include "readSIMPLEControls.H"
60         p.storePrevIter();
62         // Pressure-velocity SIMPLE corrector
63         {
64             // Momentum predictor
65             tmp<fvVectorMatrix> UEqn
66             (
67                 fvm::div(phi, U)
68               - fvm::Sp(fvc::div(phi), U)
69               + turbulence->divDevReff(U)
70             );
71             mrfZones.addCoriolis(UEqn());
73             UEqn().relax();
75             solve(UEqn() == -fvc::grad(p));
77             p.boundaryField().updateCoeffs();
78             volScalarField rAU = 1.0/UEqn().A();
79             U = rAU*UEqn().H();
80             UEqn.clear();
82             phi = fvc::interpolate(U, "interpolate(HbyA)") & mesh.Sf();
83             mrfZones.relativeFlux(phi);
84             adjustPhi(phi, U, p);
86             // Non-orthogonal pressure corrector loop
87             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
88             {
89                 fvScalarMatrix pEqn
90                 (
91                     fvm::laplacian(rAU, p) == fvc::div(phi)
92                 );
94                 pEqn.setReference(pRefCell, pRefValue);
95                 pEqn.solve();
97                 if (nonOrth == nNonOrthCorr)
98                 {
99                     phi -= pEqn.flux();
100                 }
101             }
103             #include "continuityErrs.H"
105             // Explicitly relax pressure for momentum corrector
106             p.relax();
108             // Momentum corrector
109             U -= rAU*fvc::grad(p);
110             U.correctBoundaryConditions();
111         }
113         turbulence->correct();
115         runTime.write();
117         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
118             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
119             << nl << endl;
120     }
122     Info<< "End\n" << endl;
124     return 0;
128 // ************************************************************************* //