XiFoam.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-2021 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  XiFoam
26 
27 Description
28  Solver for compressible premixed/partially-premixed combustion with
29  turbulence modelling.
30 
31  Combusting RANS code using the b-Xi two-equation model.
32  Xi may be obtained by either the solution of the Xi transport
33  equation or from an algebraic expression. Both approaches are
34  based on Gulder's flame speed correlation which has been shown
35  to be appropriate by comparison with the results from the
36  spectral model.
37 
38  Strain effects are encorporated directly into the Xi equation
39  but not in the algebraic approximation. Further work need to be
40  done on this issue, particularly regarding the enhanced removal rate
41  caused by flame compression. Analysis using results of the spectral
42  model will be required.
43 
44  For cases involving very lean Propane flames or other flames which are
45  very strain-sensitive, a transport equation for the laminar flame
46  speed is present. This equation is derived using heuristic arguments
47  involving the strain time scale and the strain-rate at extinction.
48  the transport velocity is the same as that for the Xi equation.
49 
50 \*---------------------------------------------------------------------------*/
51 
52 #include "fvCFD.H"
53 #include "psiuReactionThermo.H"
56 #include "laminarFlameSpeed.H"
57 #include "ignition.H"
58 #include "Switch.H"
59 #include "pimpleControl.H"
60 #include "fvModels.H"
61 #include "fvConstraints.H"
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 int main(int argc, char *argv[])
66 {
67  #include "postProcess.H"
68 
69  #include "setRootCaseLists.H"
70  #include "createTime.H"
71  #include "createMesh.H"
72  #include "createControl.H"
73  #include "readCombustionProperties.H"
74  #include "createFields.H"
75  #include "createFieldRefs.H"
76  #include "initContinuityErrs.H"
77  #include "createTimeControls.H"
78  #include "compressibleCourantNo.H"
79  #include "setInitialDeltaT.H"
80 
81  turbulence->validate();
82 
83  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85  Info<< "\nStarting time loop\n" << endl;
86 
87  while (pimple.run(runTime))
88  {
89  #include "readTimeControls.H"
90  #include "compressibleCourantNo.H"
91  #include "setDeltaT.H"
92 
93  runTime++;
94  Info<< "Time = " << runTime.timeName() << nl << endl;
95 
96  #include "rhoEqn.H"
97 
98  // --- Pressure-velocity PIMPLE corrector loop
99  while (pimple.loop())
100  {
101  fvModels.correct();
102 
103  #include "UEqn.H"
104 
105  #include "ftEqn.H"
106  #include "bEqn.H"
107  #include "EauEqn.H"
108  #include "EaEqn.H"
109 
110  if (!ign.ignited())
111  {
112  thermo.heu() == thermo.he();
113  }
114 
115  // --- Pressure corrector loop
116  while (pimple.correct())
117  {
118  #include "pEqn.H"
119  }
120 
121  if (pimple.turbCorr())
122  {
123  turbulence->correct();
124  thermophysicalTransport->correct();
125  }
126  }
127 
128  rho = thermo.rho();
129 
130  runTime.write();
131 
132  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
133  << " ClockTime = " << runTime.elapsedClockTime() << " s"
134  << nl << endl;
135  }
136 
137  Info<< "End\n" << endl;
138 
139  return 0;
140 }
141 
142 
143 // ************************************************************************* //
fluidReactionThermo & thermo
Definition: createFields.H:28
pimpleNoLoopControl & pimple
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:316
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Read the control parameters used by setDeltaT.
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
fluidReactionThermophysicalTransportModel & thermophysicalTransport
Foam::fvModels & fvModels
messageStream Info
Execute application functionObjects to post-process existing results.
Read the control parameters used by setDeltaT.