pisoFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Application
25  pisoFoam
26 
27 Description
28  Transient solver for incompressible flow.
29 
30  Sub-models include:
31  - turbulence modelling, i.e. laminar, RAS or LES
32  - run-time selectable MRF and finite volume options, e.g. explicit porosity
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "fvCFD.H"
39 #include "pisoControl.H"
40 #include "fvIOoptionList.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 int main(int argc, char *argv[])
45 {
46  #include "setRootCase.H"
47  #include "createTime.H"
48  #include "createMesh.H"
49 
50  pisoControl piso(mesh);
51 
52  #include "createFields.H"
53  #include "createMRF.H"
54  #include "createFvOptions.H"
55  #include "initContinuityErrs.H"
56 
57  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59  Info<< "\nStarting time loop\n" << endl;
60 
61  while (runTime.loop())
62  {
63  Info<< "Time = " << runTime.timeName() << nl << endl;
64 
65  #include "CourantNo.H"
66 
67  // Pressure-velocity PISO corrector
68  {
69  #include "UEqn.H"
70 
71  // --- PISO loop
72  while (piso.correct())
73  {
74  #include "pEqn.H"
75  }
76  }
77 
78  laminarTransport.correct();
79  turbulence->correct();
80 
81  runTime.write();
82 
83  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
84  << " ClockTime = " << runTime.elapsedClockTime() << " s"
85  << nl << endl;
86  }
87 
88  Info<< "End\n" << endl;
89 
90  return 0;
91 }
92 
93 
94 // ************************************************************************* //
singlePhaseTransportModel laminarTransport(U, phi)
messageStream Info
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
dynamicFvMesh & mesh
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
int main(int argc, char *argv[])
Definition: postCalc.C:54