ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / applications / solvers / incompressible / icoFoam / icoFoam.C
blobbfc39e572c0a6cad0f6a930c50b0d4ec4e4666df
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Application
25     icoFoam
27 Description
28     Transient solver for incompressible, laminar flow of Newtonian fluids.
30 \*---------------------------------------------------------------------------*/
32 #include "fvCFD.H"
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 int main(int argc, char *argv[])
38     #include "setRootCase.H"
40     #include "createTime.H"
41     #include "createMesh.H"
42     #include "createFields.H"
43     #include "initContinuityErrs.H"
45     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47     Info<< "\nStarting time loop\n" << endl;
49     while (runTime.loop())
50     {
51         Info<< "Time = " << runTime.timeName() << nl << endl;
53         #include "readPISOControls.H"
54         #include "CourantNo.H"
56         fvVectorMatrix UEqn
57         (
58             fvm::ddt(U)
59           + fvm::div(phi, U)
60           - fvm::laplacian(nu, U)
61         );
63         solve(UEqn == -fvc::grad(p));
65         // --- PISO loop
67         for (int corr=0; corr<nCorr; corr++)
68         {
69             volScalarField rAU(1.0/UEqn.A());
71             U = rAU*UEqn.H();
72             phi = (fvc::interpolate(U) & mesh.Sf())
73                 + fvc::ddtPhiCorr(rAU, U, phi);
75             adjustPhi(phi, U, p);
77             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
78             {
79                 fvScalarMatrix pEqn
80                 (
81                     fvm::laplacian(rAU, p) == fvc::div(phi)
82                 );
84                 pEqn.setReference(pRefCell, pRefValue);
85                 pEqn.solve();
87                 if (nonOrth == nNonOrthCorr)
88                 {
89                     phi -= pEqn.flux();
90                 }
91             }
93             #include "continuityErrs.H"
95             U -= rAU*fvc::grad(p);
96             U.correctBoundaryConditions();
97         }
99         runTime.write();
101         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
102             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
103             << nl << endl;
104     }
106     Info<< "End\n" << endl;
108     return 0;
112 // ************************************************************************* //