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