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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  XiFoam
28 
29 Group
30  grpCombustionSolvers
31 
32 Description
33  Solver for compressible premixed/partially-premixed combustion with
34  turbulence modelling.
35 
36  Combusting RANS code using the b-Xi two-equation model.
37  Xi may be obtained by either the solution of the Xi transport
38  equation or from an algebraic expression. Both approaches are
39  based on Gulder's flame speed correlation which has been shown
40  to be appropriate by comparison with the results from the
41  spectral model.
42 
43  Strain effects are encorporated directly into the Xi equation
44  but not in the algebraic approximation. Further work need to be
45  done on this issue, particularly regarding the enhanced removal rate
46  caused by flame compression. Analysis using results of the spectral
47  model will be required.
48 
49  For cases involving very lean Propane flames or other flames which are
50  very strain-sensitive, a transport equation for the laminar flame
51  speed is present. This equation is derived using heuristic arguments
52  involving the strain time scale and the strain-rate at extinction.
53  the transport velocity is the same as that for the Xi equation.
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #include "fvCFD.H"
58 #include "psiuReactionThermo.H"
60 #include "laminarFlameSpeed.H"
61 #include "ignition.H"
62 #include "Switch.H"
63 #include "pimpleControl.H"
64 #include "fvOptions.H"
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 int main(int argc, char *argv[])
69 {
70  argList::addNote
71  (
72  "Solver for compressible premixed/partially-premixed combustion with"
73  " turbulence modelling."
74  );
75 
76  #include "postProcess.H"
77 
78  #include "addCheckCaseOptions.H"
79  #include "setRootCaseLists.H"
80  #include "createTime.H"
81  #include "createMesh.H"
82  #include "createControl.H"
83  #include "readCombustionProperties.H"
84  #include "createFields.H"
85  #include "createFieldRefs.H"
86  #include "initContinuityErrs.H"
87  #include "createTimeControls.H"
88  #include "compressibleCourantNo.H"
89  #include "setInitialDeltaT.H"
90 
91  turbulence->validate();
92 
93  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
94 
95  Info<< "\nStarting time loop\n" << endl;
96 
97  while (runTime.run())
98  {
99  #include "readTimeControls.H"
100  #include "compressibleCourantNo.H"
101  #include "setDeltaT.H"
102 
103  ++runTime;
104  Info<< "Time = " << runTime.timeName() << nl << endl;
105 
106  #include "rhoEqn.H"
107 
108  // --- Pressure-velocity PIMPLE corrector loop
109  while (pimple.loop())
110  {
111  #include "UEqn.H"
112 
113  #include "ftEqn.H"
114  #include "bEqn.H"
115  #include "EauEqn.H"
116  #include "EaEqn.H"
117 
118  if (!ign.ignited())
119  {
120  thermo.heu() == thermo.he();
121  }
122 
123  // --- Pressure corrector loop
124  while (pimple.correct())
125  {
126  #include "pEqn.H"
127  }
128 
129  if (pimple.turbCorr())
130  {
131  turbulence->correct();
132  }
133  }
134 
135  rho = thermo.rho();
136 
137  runTime.write();
138 
139  runTime.printExecutionTime(Info);
140  }
141 
142  Info<< "End\n" << endl;
143 
144  return 0;
145 }
146 
147 
148 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
psiuReactionThermo.H
fvOptions.H
turbulence
Info<< "Reading field U\n"<< 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\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
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:369
laminarFlameSpeed.H
rho
rho
Definition: readInitialConditions.H:88
pimpleControl.H
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
addCheckCaseOptions.H
Required Classes.
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
ignition.H
Switch.H
postProcess.H
Execute application functionObjects to post-process existing results.
Foam::nl
constexpr char nl
Definition: Ostream.H:404
createTimeControls.H
Read the control parameters used by setDeltaT.
readTimeControls.H
Read the control parameters used by setDeltaT.
createMesh.H
Required Variables.
createTime.H
fvCFD.H
turbulentFluidThermoModel.H