ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / applications / solvers / incompressible / pisoFoam / pisoFoam.C
blob9dec6c12e90368bfdd6610c7b40f55e9a875ddae
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     pisoFoam
27 Description
28     Transient solver for incompressible flow.
30     Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
32 \*---------------------------------------------------------------------------*/
34 #include "fvCFD.H"
35 #include "singlePhaseTransportModel.H"
36 #include "turbulenceModel.H"
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 int main(int argc, char *argv[])
42     #include "setRootCase.H"
44     #include "createTime.H"
45     #include "createMesh.H"
46     #include "createFields.H"
47     #include "initContinuityErrs.H"
49     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51     Info<< "\nStarting time loop\n" << endl;
53     while (runTime.loop())
54     {
55         Info<< "Time = " << runTime.timeName() << nl << endl;
57         #include "readPISOControls.H"
58         #include "CourantNo.H"
60         // Pressure-velocity PISO corrector
61         {
62             // Momentum predictor
64             fvVectorMatrix UEqn
65             (
66                 fvm::ddt(U)
67               + fvm::div(phi, U)
68               + turbulence->divDevReff(U)
69             );
71             UEqn.relax();
73             if (momentumPredictor)
74             {
75                 solve(UEqn == -fvc::grad(p));
76             }
78             // --- PISO loop
80             for (int corr=0; corr<nCorr; corr++)
81             {
82                 volScalarField rAU(1.0/UEqn.A());
84                 U = rAU*UEqn.H();
85                 phi = (fvc::interpolate(U) & mesh.Sf())
86                     + fvc::ddtPhiCorr(rAU, U, phi);
88                 adjustPhi(phi, U, p);
90                 // Non-orthogonal pressure corrector loop
91                 for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
92                 {
93                     // Pressure corrector
95                     fvScalarMatrix pEqn
96                     (
97                         fvm::laplacian(rAU, p) == fvc::div(phi)
98                     );
100                     pEqn.setReference(pRefCell, pRefValue);
102                     if
103                     (
104                         corr == nCorr-1
105                      && nonOrth == nNonOrthCorr
106                     )
107                     {
108                         pEqn.solve(mesh.solver("pFinal"));
109                     }
110                     else
111                     {
112                         pEqn.solve();
113                     }
115                     if (nonOrth == nNonOrthCorr)
116                     {
117                         phi -= pEqn.flux();
118                     }
119                 }
121                 #include "continuityErrs.H"
123                 U -= rAU*fvc::grad(p);
124                 U.correctBoundaryConditions();
125             }
126         }
128         turbulence->correct();
130         runTime.write();
132         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
133             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
134             << nl << endl;
135     }
137     Info<< "End\n" << endl;
139     return 0;
143 // ************************************************************************* //