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
compressible
sonicFoam
sonicFoam.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) 2011-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
sonicFoam
26
27
Group
28
grpCompressibleSolvers
29
30
Description
31
Transient solver for trans-sonic/supersonic, turbulent flow of a
32
compressible gas.
33
34
\*---------------------------------------------------------------------------*/
35
36
#include "
fvCFD.H
"
37
#include "
psiThermo.H
"
38
#include "
turbulentFluidThermoModel.H
"
39
#include "
pimpleControl.H
"
40
#include "
fvOptions.H
"
41
42
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44
int
main(
int
argc,
char
*argv[])
45
{
46
#include "
postProcess.H
"
47
48
#include "
setRootCase.H
"
49
#include "
createTime.H
"
50
#include "
createMesh.H
"
51
#include "createControl.H"
52
#include "createFields.H"
53
#include "createFieldRefs.H"
54
#include "
createFvOptions.H
"
55
#include "initContinuityErrs.H"
56
57
turbulence
->validate();
58
59
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60
61
Info
<<
"\nStarting time loop\n"
<<
endl
;
62
63
while
(runTime.loop())
64
{
65
Info
<<
"Time = "
<< runTime.timeName() <<
nl
<<
endl
;
66
67
#include "compressibleCourantNo.H"
68
69
#include "rhoEqn.H"
70
71
// --- Pressure-velocity PIMPLE corrector loop
72
while
(
pimple
.loop())
73
{
74
#include "UEqn.H"
75
#include "EEqn.H"
76
77
// --- Pressure corrector loop
78
while
(
pimple
.correct())
79
{
80
#include "pEqn.H"
81
}
82
83
if
(
pimple
.turbCorr())
84
{
85
turbulence
->correct();
86
}
87
}
88
89
rho
=
thermo
.rho();
90
91
runTime.write();
92
93
Info
<<
"ExecutionTime = "
<< runTime.elapsedCpuTime() <<
" s"
94
<<
" ClockTime = "
<< runTime.elapsedClockTime() <<
" s"
95
<<
nl
<<
endl
;
96
}
97
98
Info
<<
"End\n"
<<
endl
;
99
100
return
0;
101
}
102
103
104
// ************************************************************************* //
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition:
createFields.H:23
fvCFD.H
createFvOptions.H
psiThermo.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
turbulentFluidThermoModel.H
thermo
psiReactionThermo & thermo
Definition:
createFields.H:31
rho
rho
Definition:
readInitialConditions.H:96
Foam::nl
static const char nl
Definition:
Ostream.H:262
createMesh.H
fvOptions.H
setRootCase.H
pimpleControl.H
Foam::Info
messageStream Info
postProcess.H
Execute application functionObjects to post-process existing results.
Generated by
1.8.11