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