thermoFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2018 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  thermoFoam
26 
27 Description
28  Solver for energy transport and thermodynamics on a frozen flow field.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
33 #include "rhoThermo.H"
35 #include "LESModel.H"
36 #include "radiationModel.H"
37 #include "fvOptions.H"
38 #include "simpleControl.H"
39 #include "pimpleControl.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 int main(int argc, char *argv[])
44 {
45  #define NO_CONTROL
46  #include "postProcess.H"
47 
48  #include "setRootCaseLists.H"
49  #include "createTime.H"
50  #include "createMesh.H"
51  #include "createFields.H"
52 
53  const volScalarField& alphaEff = talphaEff();
54 
55  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57  Info<< "\nEvolving thermodynamics\n" << endl;
58 
59  if (mesh.solutionDict().found("SIMPLE"))
60  {
61  simpleControl simple(mesh);
62 
63  while (simple.loop(runTime))
64  {
65  Info<< "Time = " << runTime.timeName() << nl << endl;
66 
67  while (simple.correctNonOrthogonal())
68  {
69  #include "EEqn.H"
70  }
71 
72  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
73  << " ClockTime = " << runTime.elapsedClockTime() << " s"
74  << nl << endl;
75 
76  runTime.write();
77  }
78  }
79  else
80  {
81  pimpleControl pimple(mesh);
82 
83  while (runTime.run())
84  {
85  runTime++;
86 
87  Info<< "Time = " << runTime.timeName() << nl << endl;
88 
89  while (pimple.correctNonOrthogonal())
90  {
91  #include "EEqn.H"
92  }
93 
94  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
95  << " ClockTime = " << runTime.elapsedClockTime() << " s"
96  << nl << endl;
97 
98  runTime.write();
99  }
100  }
101 
102  Info<< "End\n" << endl;
103 
104  return 0;
105 }
106 
107 
108 // ************************************************************************* //
pimpleNoLoopControl & pimple
Info<< "Creating turbulence model\"<< endl;tmp< volScalarField > talphaEff
Definition: setAlphaEff.H:2
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
static const char nl
Definition: Ostream.H:265
messageStream Info
Execute application functionObjects to post-process existing results.
simpleControl simple(mesh)