OpenFOAM
8
The OpenFOAM Foundation
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
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-2020 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
"
34
#include "
fluidThermoMomentumTransportModel.H
"
35
#include "
fluidThermophysicalTransportModel.H
"
36
#include "
LESModel.H
"
37
#include "
radiationModel.H
"
38
#include "
fvOptions.H
"
39
#include "
simpleControl.H
"
40
#include "
pimpleControl.H
"
41
42
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44
int
main(
int
argc,
char
*argv[])
45
{
46
#define NO_CONTROL
47
#include "
postProcess.H
"
48
49
#include "
setRootCaseLists.H
"
50
#include "
createTime.H
"
51
#include "
createMesh.H
"
52
#include "createFields.H"
53
54
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55
56
Info
<<
"\nEvolving thermodynamics\n"
<<
endl
;
57
58
if
(
mesh
.solutionDict().found(
"SIMPLE"
))
59
{
60
simpleControl
simple
(
mesh
);
61
62
while
(
simple
.loop(
runTime
))
63
{
64
Info
<<
"Time = "
<<
runTime
.timeName() <<
nl
<<
endl
;
65
66
while
(
simple
.correctNonOrthogonal())
67
{
68
#include "EEqn.H"
69
}
70
71
Info
<<
"ExecutionTime = "
<<
runTime
.elapsedCpuTime() <<
" s"
72
<<
" ClockTime = "
<<
runTime
.elapsedClockTime() <<
" s"
73
<<
nl
<<
endl
;
74
75
runTime
.write();
76
}
77
}
78
else
79
{
80
pimpleControl
pimple
(
mesh
);
81
82
while
(
runTime
.run())
83
{
84
runTime
++;
85
86
Info
<<
"Time = "
<<
runTime
.timeName() <<
nl
<<
endl
;
87
88
while
(
pimple
.correctNonOrthogonal())
89
{
90
#include "EEqn.H"
91
}
92
93
Info
<<
"ExecutionTime = "
<<
runTime
.elapsedCpuTime() <<
" s"
94
<<
" ClockTime = "
<<
runTime
.elapsedClockTime() <<
" s"
95
<<
nl
<<
endl
;
96
97
runTime
.write();
98
}
99
}
100
101
Info
<<
"End\n"
<<
endl
;
102
103
return
0;
104
}
105
106
107
// ************************************************************************* //
fvCFD.H
pimple
pimpleNoLoopControl & pimple
Definition:
setRegionFluidFields.H:62
setRootCaseLists.H
radiationModel.H
runTime
engineTime & runTime
Definition:
createEngineTime.H:13
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition:
Ostream.H:251
createTime.H
LESModel.H
mesh
dynamicFvMesh & mesh
Definition:
createDynamicFvMesh.H:18
Foam::nl
static const char nl
Definition:
Ostream.H:260
fluidThermophysicalTransportModel.H
createMesh.H
fvOptions.H
pimpleControl.H
Foam::Info
messageStream Info
postProcess.H
Execute application functionObjects to post-process existing results.
rhoThermo.H
fluidThermoMomentumTransportModel.H
simpleControl.H
simple
simpleControl simple(mesh)
applications
solvers
heatTransfer
thermoFoam
thermoFoam.C
Generated by
1.8.13