initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / solvers / incompressible / icoFoam / icoFoam.C
blobfef690f831a08a6db3da3c936caa96c35ac5c445
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     icoFoam
28 Description
29     Transient solver for incompressible, laminar flow of Newtonian fluids.
31 \*---------------------------------------------------------------------------*/
33 #include "fvCFD.H"
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 int main(int argc, char *argv[])
39     #include "setRootCase.H"
41     #include "createTime.H"
42     #include "createMesh.H"
43     #include "createFields.H"
44     #include "initContinuityErrs.H"
46     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48     Info<< "\nStarting time loop\n" << endl;
50     while (runTime.loop())
51     {
52         Info<< "Time = " << runTime.timeName() << nl << endl;
54         #include "readPISOControls.H"
55         #include "CourantNo.H"
57         fvVectorMatrix UEqn
58         (
59             fvm::ddt(U)
60           + fvm::div(phi, U)
61           - fvm::laplacian(nu, U)
62         );
64         solve(UEqn == -fvc::grad(p));
66         // --- PISO loop
68         for (int corr=0; corr<nCorr; corr++)
69         {
70             volScalarField rUA = 1.0/UEqn.A();
72             U = rUA*UEqn.H();
73             phi = (fvc::interpolate(U) & mesh.Sf())
74                 + fvc::ddtPhiCorr(rUA, U, phi);
76             adjustPhi(phi, U, p);
78             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
79             {
80                 fvScalarMatrix pEqn
81                 (
82                     fvm::laplacian(rUA, p) == fvc::div(phi)
83                 );
85                 pEqn.setReference(pRefCell, pRefValue);
86                 pEqn.solve();
88                 if (nonOrth == nNonOrthCorr)
89                 {
90                     phi -= pEqn.flux();
91                 }
92             }
94             #include "continuityErrs.H"
96             U -= rUA*fvc::grad(p);
97             U.correctBoundaryConditions();
98         }
100         runTime.write();
102         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
103             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
104             << nl << endl;
105     }
107     Info<< "End\n" << endl;
109     return 0;
113 // ************************************************************************* //