simpleFoam.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  simpleFoam
28 
29 Group
30  grpIncompressibleSolvers
31 
32 Description
33  Steady-state solver for incompressible, turbulent flows.
34 
35  \heading Solver details
36  The solver uses the SIMPLE algorithm to solve the continuity equation:
37 
38  \f[
39  \div \vec{U} = 0
40  \f]
41 
42  and momentum equation:
43 
44  \f[
45  \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  \heading Required fields
58  \plaintable
59  U | Velocity [m/s]
60  p | Kinematic pressure, p/rho [m2/s2]
61  <turbulence fields> | As required by user selection
62  \endplaintable
63 
64 \*---------------------------------------------------------------------------*/
65 
66 #include "fvCFD.H"
67 #include "dynamicFvMesh.H"
70 #include "simpleControl.H"
71 #include "fvOptions.H"
72 
73 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75 int main(int argc, char *argv[])
76 {
77  argList::addNote
78  (
79  "Steady-state solver for incompressible, turbulent flows."
80  );
81 
82  #include "postProcess.H"
83 
84  #include "addCheckCaseOptions.H"
85  #include "setRootCaseLists.H"
86  #include "createTime.H"
87  #include "createDynamicFvMesh.H"
88  #include "createControl.H"
89  #include "createFields.H"
90  #include "initContinuityErrs.H"
91 
92  turbulence->validate();
93 
94  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
95 
96  Info<< "\nStarting time loop\n" << endl;
97 
98  while (simple.loop())
99  {
100  Info<< "Time = " << runTime.timeName() << nl << endl;
101 
102  // Do any mesh changes
103  mesh.controlledUpdate();
104 
105  if (mesh.changing())
106  {
107  MRF.update();
108  }
109 
110  // --- Pressure-velocity SIMPLE corrector
111  {
112  #include "UEqn.H"
113  #include "pEqn.H"
114  }
115 
116  laminarTransport.correct();
117  turbulence->correct();
118 
119  runTime.write();
120 
121  runTime.printExecutionTime(Info);
122  }
123 
124  Info<< "End\n" << endl;
125 
126  return 0;
127 }
128 
129 
130 // ************************************************************************* //
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
singlePhaseTransportModel laminarTransport(U, phi)
dynamicFvMesh & mesh
IOMRFZoneList & MRF
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
const dictionary & simple
messageStream Info
Information stream (stdout output on master, null elsewhere)
Execute application functionObjects to post-process existing results.
Required Classes.