OpenFOAM
4.1
The OpenFOAM Foundation
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Modules
Pages
applications
solvers
heatTransfer
thermoFoam
thermoFoam.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) 2013-2016 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 "
turbulentFluidThermoModel.H
"
35
#include "
turbulentFluidThermoModel.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 "
setRootCase.H
"
50
#include "
createTime.H
"
51
#include "
createMesh.H
"
52
#include "createFields.H"
53
#include "
createFvOptions.H
"
54
55
const
volScalarField
& alphaEff =
talphaEff
();
56
57
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58
59
Info
<<
"\nEvolving thermodynamics\n"
<<
endl
;
60
61
if
(
mesh
.solutionDict().found(
"SIMPLE"
))
62
{
63
simpleControl
simple
(
mesh
);
64
65
while
(
simple
.loop())
66
{
67
Info
<<
"Time = "
<< runTime.timeName() <<
nl
<<
endl
;
68
69
while
(
simple
.correctNonOrthogonal())
70
{
71
#include "EEqn.H"
72
}
73
74
Info
<<
"ExecutionTime = "
<< runTime.elapsedCpuTime() <<
" s"
75
<<
" ClockTime = "
<< runTime.elapsedClockTime() <<
" s"
76
<<
nl
<<
endl
;
77
78
runTime.
write
();
79
}
80
}
81
else
82
{
83
pimpleControl
pimple
(
mesh
);
84
85
while
(runTime.run())
86
{
87
runTime++;
88
89
Info
<<
"Time = "
<< runTime.timeName() <<
nl
<<
endl
;
90
91
while
(
pimple
.correctNonOrthogonal())
92
{
93
#include "EEqn.H"
94
}
95
96
Info
<<
"ExecutionTime = "
<< runTime.elapsedCpuTime() <<
" s"
97
<<
" ClockTime = "
<< runTime.elapsedClockTime() <<
" s"
98
<<
nl
<<
endl
;
99
100
runTime.
write
();
101
}
102
}
103
104
Info
<<
"End\n"
<<
endl
;
105
106
return
0;
107
}
108
109
110
// ************************************************************************* //
fvCFD.H
createFvOptions.H
talphaEff
Info<< "Creating turbulence model\n"<< endl;tmp< volScalarField > talphaEff
Definition:
setAlphaEff.H:2
radiationModel.H
pimple
const dictionary & pimple
Definition:
readFluidMultiRegionPIMPLEControls.H:1
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition:
Ostream.H:253
createTime.H
LESModel.H
turbulentFluidThermoModel.H
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition:
volFieldsFwd.H:52
mesh
dynamicFvMesh & mesh
Definition:
createDynamicFvMesh.H:18
Foam::nl
static const char nl
Definition:
Ostream.H:262
createMesh.H
fvOptions.H
setRootCase.H
simple
const dictionary & simple
Definition:
readFluidMultiRegionSIMPLEControls.H:1
pimpleControl.H
Foam::Info
messageStream Info
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
postProcess.H
Execute application functionObjects to post-process existing results.
rhoThermo.H
simpleControl.H
Generated by
1.8.11