reactingFoam.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  reactingFoam
28 
29 Group
30  grpCombustionSolvers
31 
32 Description
33  Solver for combustion with chemical reactions.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "fvCFD.H"
39 #include "psiReactionThermo.H"
40 #include "CombustionModel.H"
41 #include "multivariateScheme.H"
42 #include "pimpleControl.H"
43 #include "pressureControl.H"
44 #include "fvOptions.H"
45 #include "localEulerDdtScheme.H"
46 #include "fvcSmooth.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 int main(int argc, char *argv[])
51 {
52  argList::addNote
53  (
54  "Solver for combustion with chemical reactions"
55  );
56 
57  #include "postProcess.H"
58 
59  #include "addCheckCaseOptions.H"
60  #include "setRootCaseLists.H"
61  #include "createTime.H"
62  #include "createMesh.H"
63  #include "createControl.H"
64  #include "createTimeControls.H"
65  #include "initContinuityErrs.H"
66  #include "createFields.H"
67  #include "createFieldRefs.H"
68 
69  turbulence->validate();
70 
71  if (!LTS)
72  {
73  #include "compressibleCourantNo.H"
74  #include "setInitialDeltaT.H"
75  }
76 
77  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78 
79  Info<< "\nStarting time loop\n" << endl;
80 
81  while (runTime.run())
82  {
83  #include "readTimeControls.H"
84 
85  if (LTS)
86  {
87  #include "setRDeltaT.H"
88  }
89  else
90  {
91  #include "compressibleCourantNo.H"
92  #include "setDeltaT.H"
93  }
94 
95  ++runTime;
96 
97  Info<< "Time = " << runTime.timeName() << nl << endl;
98 
99  #include "rhoEqn.H"
100 
101  while (pimple.loop())
102  {
103  #include "UEqn.H"
104  #include "YEqn.H"
105  #include "EEqn.H"
106 
107  // --- Pressure corrector loop
108  while (pimple.correct())
109  {
110  if (pimple.consistent())
111  {
112  #include "pcEqn.H"
113  }
114  else
115  {
116  #include "pEqn.H"
117  }
118  }
119 
120  if (pimple.turbCorr())
121  {
122  turbulence->correct();
123  }
124  }
125 
126  rho = thermo.rho();
127 
128  runTime.write();
129 
130  runTime.printExecutionTime(Info);
131  }
132 
133  Info<< "End\n" << endl;
134 
135  return 0;
136 }
137 
138 
139 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
fvOptions.H
multivariateScheme.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
rho
rho
Definition: readInitialConditions.H:88
pimpleControl.H
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
addCheckCaseOptions.H
Required Classes.
localEulerDdtScheme.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
CombustionModel.H
LTS
bool LTS
Definition: createRDeltaT.H:1
postProcess.H
Execute application functionObjects to post-process existing results.
pressureControl.H
psiReactionThermo.H
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
fvcSmooth.H
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...