engineFoam.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  engineFoam
26 
27 Description
28  Solver for internal combustion engines.
29 
30  Combusting RANS code using the b-Xi two-equation model.
31  Xi may be obtained by either the solution of the Xi transport
32  equation or from an algebraic exression. Both approaches are
33  based on Gulder's flame speed correlation which has been shown
34  to be appropriate by comparison with the results from the
35  spectral model.
36 
37  Strain effects are encorporated directly into the Xi equation
38  but not in the algebraic approximation. Further work need to be
39  done on this issue, particularly regarding the enhanced removal rate
40  caused by flame compression. Analysis using results of the spectral
41  model will be required.
42 
43  For cases involving very lean Propane flames or other flames which are
44  very strain-sensitive, a transport equation for the laminar flame
45  speed is present. This equation is derived using heuristic arguments
46  involving the strain time scale and the strain-rate at extinction.
47  the transport velocity is the same as that for the Xi equation.
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #include "fvCFD.H"
52 #include "engineTime.H"
53 #include "engineMesh.H"
54 #include "psiuReactionThermo.H"
56 #include "laminarFlameSpeed.H"
57 #include "ignition.H"
58 #include "Switch.H"
59 #include "OFstream.H"
60 #include "mathematicalConstants.H"
61 #include "pimpleControl.H"
62 #include "fvIOoptionList.H"
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 int main(int argc, char *argv[])
67 {
68  #include "setRootCase.H"
69 
70  #include "createEngineTime.H"
71  #include "createEngineMesh.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 "createRhoUf.H"
80  #include "initContinuityErrs.H"
81  #include "readEngineTimeControls.H"
82  #include "compressibleCourantNo.H"
83  #include "setInitialDeltaT.H"
84  #include "startSummary.H"
85 
86  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
87 
88  Info<< "\nStarting time loop\n" << endl;
89 
90  while (runTime.run())
91  {
92  #include "readEngineTimeControls.H"
93  #include "compressibleCourantNo.H"
94  #include "setDeltaT.H"
95 
96  runTime++;
97 
98  Info<< "Crank angle = " << runTime.theta() << " CA-deg" << endl;
99 
100  mesh.move();
101 
102  #include "rhoEqn.H"
103 
104  // --- Pressure-velocity PIMPLE corrector loop
105  while (pimple.loop())
106  {
107  #include "UEqn.H"
108 
109  #include "ftEqn.H"
110  #include "bEqn.H"
111  #include "EauEqn.H"
112  #include "EaEqn.H"
113 
114  if (!ign.ignited())
115  {
116  thermo.heu() == thermo.he();
117  }
118 
119  // --- Pressure corrector loop
120  while (pimple.correct())
121  {
122  #include "pEqn.H"
123  }
124 
125  if (pimple.turbCorr())
126  {
127  turbulence->correct();
128  }
129  }
130 
131  #include "logSummary.H"
132 
133  rho = thermo.rho();
134 
135  runTime.write();
136 
137  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
138  << " ClockTime = " << runTime.elapsedClockTime() << " s"
139  << nl << endl;
140  }
141 
142  Info<< "End\n" << endl;
143 
144  return 0;
145 }
146 
147 
148 // ************************************************************************* //
messageStream Info
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
dynamicFvMesh & mesh
const dictionary & pimple
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
int main(int argc, char *argv[])
Definition: postCalc.C:54
psiReactionThermo & thermo
Definition: createFields.H:32
Creates and initialises the velocity velocity field rhoUf.