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
combustion
XiFoam
XiFoam.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
XiFoam
26
27
Description
28
Solver for compressible premixed/partially-premixed combustion with
29
turbulence modelling.
30
31
Combusting RANS code using the b-Xi two-equation model.
32
Xi may be obtained by either the solution of the Xi transport
33
equation or from an algebraic exression. Both approaches are
34
based on Gulder's flame speed correlation which has been shown
35
to be appropriate by comparison with the results from the
36
spectral model.
37
38
Strain effects are encorporated directly into the Xi equation
39
but not in the algebraic approximation. Further work need to be
40
done on this issue, particularly regarding the enhanced removal rate
41
caused by flame compression. Analysis using results of the spectral
42
model will be required.
43
44
For cases involving very lean Propane flames or other flames which are
45
very strain-sensitive, a transport equation for the laminar flame
46
speed is present. This equation is derived using heuristic arguments
47
involving the strain time scale and the strain-rate at extinction.
48
the transport velocity is the same as that for the Xi equation.
49
50
\*---------------------------------------------------------------------------*/
51
52
#include "
fvCFD.H
"
53
#include "
psiuReactionThermo.H
"
54
#include "
turbulentFluidThermoModel.H
"
55
#include "
laminarFlameSpeed.H
"
56
#include "
ignition.H
"
57
#include "
Switch.H
"
58
#include "
pimpleControl.H
"
59
#include "
fvOptions.H
"
60
61
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62
63
int
main(
int
argc,
char
*argv[])
64
{
65
#include "
postProcess.H
"
66
67
#include "
setRootCase.H
"
68
#include "
createTime.H
"
69
#include "
createMesh.H
"
70
#include "createControl.H"
71
#include "readCombustionProperties.H"
72
#include "createFields.H"
73
#include "createFieldRefs.H"
74
#include "
createFvOptions.H
"
75
#include "initContinuityErrs.H"
76
#include "
createTimeControls.H
"
77
#include "compressibleCourantNo.H"
78
#include "setInitialDeltaT.H"
79
80
turbulence
->validate();
81
82
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83
84
Info
<<
"\nStarting time loop\n"
<<
endl
;
85
86
while
(runTime.run())
87
{
88
#include "
readTimeControls.H
"
89
#include "compressibleCourantNo.H"
90
#include "setDeltaT.H"
91
92
runTime++;
93
Info
<<
"Time = "
<< runTime.timeName() <<
nl
<<
endl
;
94
95
#include "rhoEqn.H"
96
97
// --- Pressure-velocity PIMPLE corrector loop
98
while
(
pimple
.loop())
99
{
100
#include "UEqn.H"
101
102
#include "ftEqn.H"
103
#include "bEqn.H"
104
#include "EauEqn.H"
105
#include "EaEqn.H"
106
107
if
(!ign.ignited())
108
{
109
thermo
.heu() ==
thermo
.he();
110
}
111
112
// --- Pressure corrector loop
113
while
(
pimple
.correct())
114
{
115
#include "pEqn.H"
116
}
117
118
if
(
pimple
.turbCorr())
119
{
120
turbulence
->correct();
121
}
122
}
123
124
rho
=
thermo
.rho();
125
126
runTime.write();
127
128
Info
<<
"ExecutionTime = "
<< runTime.elapsedCpuTime() <<
" s"
129
<<
" ClockTime = "
<< runTime.elapsedClockTime() <<
" s"
130
<<
nl
<<
endl
;
131
}
132
133
Info
<<
"End\n"
<<
endl
;
134
135
return
0;
136
}
137
138
139
// ************************************************************************* //
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition:
createFields.H:23
fvCFD.H
createFvOptions.H
pimple
const dictionary & pimple
Definition:
readFluidMultiRegionPIMPLEControls.H:1
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition:
Ostream.H:253
readTimeControls.H
Read the control parameters used by setDeltaT.
createTime.H
laminarFlameSpeed.H
turbulentFluidThermoModel.H
thermo
psiReactionThermo & thermo
Definition:
createFields.H:31
rho
rho
Definition:
readInitialConditions.H:96
psiuReactionThermo.H
Foam::nl
static const char nl
Definition:
Ostream.H:262
createMesh.H
fvOptions.H
setRootCase.H
pimpleControl.H
Foam::Info
messageStream Info
ignition.H
Switch.H
postProcess.H
Execute application functionObjects to post-process existing results.
createTimeControls.H
Read the control parameters used by setDeltaT.
Generated by
1.8.11