overSimpleFoam.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-2014 OpenFOAM Foundation
9  Copyright (C) 2016-2017 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  overSimpleFoam
29 
30 Group
31  grpIncompressibleSolvers
32 
33 Description
34  Steady-state solver for incompressible flows with turbulence modelling.
35 
36  \heading Solver details
37  The solver uses the SIMPLE algorithm to solve the continuity equation:
38 
39  \f[
40  \div \vec{U} = 0
41  \f]
42 
43  and momentum equation:
44 
45  \f[
46  \div \left( \vec{U} \vec{U} \right) - \div \gvec{R}
47  = - \grad p + \vec{S}_U
48  \f]
49 
50  Where:
51  \vartable
52  \vec{U} | Velocity
53  p | Pressure
54  \vec{R} | Stress tensor
55  \vec{S}_U | Momentum source
56  \endvartable
57 
58  \heading Required fields
59  \plaintable
60  U | Velocity [m/s]
61  p | Kinematic pressure, p/rho [m2/s2]
62  <turbulence fields> | As required by user selection
63  \endplaintable
64 
65 \*---------------------------------------------------------------------------*/
66 
67 #include "fvCFD.H"
70 #include "simpleControl.H"
71 #include "fvOptions.H"
72 
73 #include "dynamicFvMesh.H"
74 #include "cellCellStencilObject.H"
75 #include "localMin.H"
76 #include "interpolationCellPoint.H"
77 #include "fvMeshSubset.H"
78 #include "transform.H"
79 #include "oversetAdjustPhi.H"
80 
81 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83 int main(int argc, char *argv[])
84 {
85  argList::addNote
86  (
87  "Steady-state solver for incompressible, turbulent flow"
88  );
89 
90  #define CREATE_MESH createUpdatedDynamicFvMesh.H
91  #include "postProcess.H"
92 
93  #include "setRootCaseLists.H"
94  #include "createTime.H"
95  #include "createUpdatedDynamicFvMesh.H"
96  #include "createControl.H"
97  #include "createFields.H"
98  #include "createOversetFields.H"
99  #include "createFvOptions.H"
100  #include "initContinuityErrs.H"
101 
102  turbulence->validate();
103 
104  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
105 
106  Info<< "\nStarting time loop\n" << endl;
107 
108  while (simple.loop())
109  {
110  Info<< "Time = " << runTime.timeName() << nl << endl;
111 
112  // --- Pressure-velocity SIMPLE corrector
113  {
114  #include "UEqn.H"
115  #include "pEqn.H"
116  }
117 
118  laminarTransport.correct();
119  turbulence->correct();
120 
121  runTime.write();
122 
123  runTime.printExecutionTime(Info);
124  }
125 
126  Info<< "End\n" << endl;
127 
128  return 0;
129 }
130 
131 
132 // ************************************************************************* //
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
singlePhaseTransportModel laminarTransport(U, phi)
3D tensor transformation operations.
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.
Adjust the balance of fluxes on the faces between interpolated and calculated to obey continuity...