pimpleFoam.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  pimpleFoam
26 
27 Group
28  grpIncompressibleSolvers
29 
30 Description
31  Large time-step transient solver for incompressible, flow using the PIMPLE
32  (merged PISO-SIMPLE) algorithm.
33 
34  \heading Solver details
35  The solver uses the PIMPLE (merged PISO-SIMPLE) algorithm to solve the
36  continuity equation:
37 
38  \f[
39  \div \vec{U} = 0
40  \f]
41 
42  and momentum equation:
43 
44  \f[
45  \ddt{\vec{U}} + \div \left( \vec{U} \vec{U} \right) - \div \gvec{R}
46  = - \grad p + \vec{S}_U
47  \f]
48 
49  Where:
50  \vartable
51  \vec{U} | Velocity
52  p | Pressure
53  \vec{R} | Stress tensor
54  \vec{S}_U | Momentum source
55  \endvartable
56 
57  Sub-models include:
58  - turbulence modelling, i.e. laminar, RAS or LES
59  - run-time selectable MRF and finite volume options, e.g. explicit porosity
60 
61  \heading Required fields
62  \plaintable
63  U | Velocity [m/s]
64  p | Kinematic pressure, p/rho [m2/s2]
65  <turbulence fields> | As required by user selection
66  \endplaintable
67 
68 \*---------------------------------------------------------------------------*/
69 
70 #include "fvCFD.H"
73 #include "pimpleControl.H"
74 #include "fvOptions.H"
76 
77 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78 
79 int main(int argc, char *argv[])
80 {
81  #include "setRootCase.H"
82  #include "createTime.H"
83  #include "createMesh.H"
84 
85  pimpleControl pimple(mesh);
86 
87  #include "createTimeControls.H"
88  #include "createFields.H"
89  #include "createMRF.H"
90  #include "createFvOptions.H"
91  #include "initContinuityErrs.H"
92 
93  turbulence->validate();
94 
95  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
96 
97  Info<< "\nStarting time loop\n" << endl;
98 
99  while (runTime.run())
100  {
101  #include "readTimeControls.H"
102  #include "CourantNo.H"
103  #include "setDeltaT.H"
104 
105  runTime++;
106 
107  Info<< "Time = " << runTime.timeName() << nl << endl;
108 
109  // --- Pressure-velocity PIMPLE corrector loop
110  while (pimple.loop())
111  {
112  #include "UEqn.H"
113 
114  // --- Pressure corrector loop
115  while (pimple.correct())
116  {
117  #include "pEqn.H"
118  }
119 
120  if (pimple.turbCorr())
121  {
122  laminarTransport.correct();
123  turbulence->correct();
124  }
125  }
126 
127  runTime.write();
128 
129  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
130  << " ClockTime = " << runTime.elapsedClockTime() << " s"
131  << nl << endl;
132  }
133 
134  Info<< "End\n" << endl;
135 
136  return 0;
137 }
138 
139 
140 // ************************************************************************* //
fvOptions.H
singlePhaseTransportModel.H
turbulentTransportModel.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
pimpleControl.H
createFvOptions.H
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
setRootCase.H
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
createTimeControls.H
Read the control parameters used by setDeltaT.
readTimeControls.H
Read the control parameters used by setDeltaT.
createMesh.H
createTime.H
fvCFD.H
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
fixedFluxPressureFvPatchScalarField.H