XiEngineFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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  XiEngineFoam
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 "fvOptions.H"
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 int main(int argc, char *argv[])
67 {
68  #define CREATE_TIME createEngineTime.H
69  #define CREATE_MESH createEngineMesh.H
70  #include "postProcess.H"
71 
72  #include "setRootCaseLists.H"
73  #include "createEngineTime.H"
74  #include "createEngineMesh.H"
75  #include "createControl.H"
76  #include "readCombustionProperties.H"
77  #include "createFields.H"
78  #include "createFieldRefs.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  turbulence->validate();
87 
88  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
89 
90  Info<< "\nStarting time loop\n" << endl;
91 
92  while (runTime.run())
93  {
94  #include "readEngineTimeControls.H"
95  #include "compressibleCourantNo.H"
96  #include "setDeltaT.H"
97 
98  runTime++;
99 
100  Info<< "Crank angle = " << runTime.theta() << " CA-deg" << endl;
101 
102  mesh.move();
103 
104  #include "rhoEqn.H"
105 
106  // --- Pressure-velocity PIMPLE corrector loop
107  while (pimple.loop())
108  {
109  #include "UEqn.H"
110 
111  #include "ftEqn.H"
112  #include "bEqn.H"
113  #include "EauEqn.H"
114  #include "EaEqn.H"
115 
116  if (!ign.ignited())
117  {
118  thermo.heu() == thermo.he();
119  }
120 
121  // --- Pressure corrector loop
122  while (pimple.correct())
123  {
124  #include "pEqn.H"
125  }
126 
127  if (pimple.turbCorr())
128  {
129  turbulence->correct();
130  }
131  }
132 
133  #include "logSummary.H"
134 
135  rho = thermo.rho();
136 
137  runTime.write();
138 
139  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
140  << " ClockTime = " << runTime.elapsedClockTime() << " s"
141  << nl << endl;
142  }
143 
144  Info<< "End\n" << endl;
145 
146  return 0;
147 }
148 
149 
150 // ************************************************************************* //
pimpleNoLoopControl & pimple
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
rhoReactionThermo & thermo
Definition: createFields.H:28
dynamicFvMesh & mesh
static const char nl
Definition: Ostream.H:265
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
messageStream Info
Execute application functionObjects to post-process existing results.
Creates and initialises the velocity velocity field rhoUf.