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-2022 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 "CorrectPhi.H"
61 #include "fvModels.H"
62 #include "fvConstraints.H"
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 int main(int argc, char *argv[])
67 {
68  #include "postProcess.H"
69 
70  #include "setRootCaseLists.H"
71  #include "createTime.H"
72  #include "createMesh.H"
73  #include "createDyMControls.H"
74  #include "initContinuityErrs.H"
75  #include "readCombustionProperties.H"
76  #include "createFields.H"
77  #include "createFieldRefs.H"
78  #include "createRhoUfIfPresent.H"
79 
80  turbulence->validate();
81 
82  #include "compressibleCourantNo.H"
83  #include "setInitialDeltaT.H"
84 
85  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87  Info<< "\nStarting time loop\n" << endl;
88 
89  while (pimple.run(runTime))
90  {
91  #include "readDyMControls.H"
92 
93  // Store divrhoU from the previous mesh so that it can be mapped
94  // and used in correctPhi to ensure the corrected phi has the
95  // same divergence
96  autoPtr<volScalarField> divrhoU;
97  if (correctPhi)
98  {
99  divrhoU = new volScalarField
100  (
101  "divrhoU",
103  );
104  }
105 
106  #include "compressibleCourantNo.H"
107  #include "setDeltaT.H"
108 
110 
111  // Store momentum to set rhoUf for introduced faces.
112  autoPtr<volVectorField> rhoU;
113  if (rhoUf.valid())
114  {
115  rhoU = new volVectorField("rhoU", rho*U);
116  }
117 
118  // Update the mesh for topology change, mesh to mesh mapping
119  mesh.update();
120 
121  runTime++;
122 
123  Info<< "Time = " << runTime.userTimeName() << nl << endl;
124 
125  // --- Pressure-velocity PIMPLE corrector loop
126  while (pimple.loop())
127  {
128  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
129  {
130  // Move the mesh
131  mesh.move();
132 
133  if (mesh.changing())
134  {
135  MRF.update();
136 
137  if (correctPhi)
138  {
139  #include "correctPhi.H"
140  }
141 
142  if (checkMeshCourantNo)
143  {
144  #include "meshCourantNo.H"
145  }
146  }
147  }
148 
149  if
150  (
151  !pimple.simpleRho()
152  && pimple.firstPimpleIter()
153  )
154  {
155  #include "rhoEqn.H"
156  }
157 
158  fvModels.correct();
159 
160  #include "UEqn.H"
161  #include "ftEqn.H"
162  #include "bEqn.H"
163  #include "EauEqn.H"
164  #include "EaEqn.H"
165 
166  if (!ign.ignited())
167  {
168  thermo.heu() == thermo.he();
169  }
170 
171  // --- Pressure corrector loop
172  while (pimple.correct())
173  {
174  #include "pEqn.H"
175  }
176 
177  if (pimple.turbCorr())
178  {
179  turbulence->correct();
180  thermophysicalTransport->correct();
181  }
182  }
183 
184  rho = thermo.rho();
185 
186  runTime.write();
187 
188  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
189  << " ClockTime = " << runTime.elapsedClockTime() << " s"
190  << nl << endl;
191  }
192 
193  Info<< "End\n" << endl;
194 
195  return 0;
196 }
197 
198 
199 // ************************************************************************* //
fluidReactionThermo & thermo
Definition: createFields.H:28
pimpleNoLoopControl & pimple
U
Definition: pEqn.H:72
IOMRFZoneList & MRF
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:353
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:255
correctPhi
checkMeshCourantNo
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Creates and initialises the velocity field rhoUf if required.
phi
Definition: correctPhi.H:3
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
moveMeshOuterCorrectors
Foam::fvModels & fvModels
Calculates and outputs the mean and maximum Courant Numbers.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
messageStream Info
Execute application functionObjects to post-process existing results.
autoPtr< surfaceVectorField > rhoUf