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-2020 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 expression. 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"
57 #include "laminarFlameSpeed.H"
58 #include "ignition.H"
59 #include "Switch.H"
60 #include "OFstream.H"
61 #include "mathematicalConstants.H"
62 #include "pimpleControl.H"
63 #include "fvOptions.H"
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 int main(int argc, char *argv[])
68 {
69  #define CREATE_TIME createEngineTime.H
70  #define CREATE_MESH createEngineMesh.H
71  #include "postProcess.H"
72 
73  #include "setRootCaseLists.H"
74  #include "createEngineTime.H"
75  #include "createEngineMesh.H"
76  #include "createControl.H"
77  #include "readCombustionProperties.H"
78  #include "createFields.H"
79  #include "createFieldRefs.H"
80  #include "createRhoUf.H"
81  #include "initContinuityErrs.H"
82  #include "readEngineTimeControls.H"
83  #include "compressibleCourantNo.H"
84  #include "setInitialDeltaT.H"
85  #include "startSummary.H"
86 
87  turbulence->validate();
88 
89  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
90 
91  Info<< "\nStarting time loop\n" << endl;
92 
93  while (pimple.run(runTime))
94  {
95  #include "readEngineTimeControls.H"
96  #include "compressibleCourantNo.H"
97  #include "setDeltaT.H"
98 
99  runTime++;
100 
101  Info<< "Crank angle = " << runTime.theta() << " CA-deg" << endl;
102 
103  mesh.move();
104 
105  #include "rhoEqn.H"
106 
107  // --- Pressure-velocity PIMPLE corrector loop
108  while (pimple.loop())
109  {
110  #include "UEqn.H"
111 
112  #include "ftEqn.H"
113  #include "bEqn.H"
114  #include "EauEqn.H"
115  #include "EaEqn.H"
116 
117  if (!ign.ignited())
118  {
119  thermo.heu() == thermo.he();
120  }
121 
122  // --- Pressure corrector loop
123  while (pimple.correct())
124  {
125  #include "pEqn.H"
126  }
127 
128  if (pimple.turbCorr())
129  {
130  turbulence->correct();
131  thermophysicalTransport->correct();
132  }
133  }
134 
135  #include "logSummary.H"
136 
137  rho = thermo.rho();
138 
139  runTime.write();
140 
141  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
142  << " ClockTime = " << runTime.elapsedClockTime() << " s"
143  << nl << endl;
144  }
145 
146  Info<< "End\n" << endl;
147 
148  return 0;
149 }
150 
151 
152 // ************************************************************************* //
pimpleNoLoopControl & pimple
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
rhoReactionThermo & thermo
Definition: createFields.H:28
dynamicFvMesh & mesh
rhoReactionThermophysicalTransportModel & thermophysicalTransport
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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
static const char nl
Definition: Ostream.H:260
messageStream Info
Execute application functionObjects to post-process existing results.
Creates and initialises the velocity velocity field rhoUf.