All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
reactingFoam.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  reactingFoam
26 
27 Description
28  Solver for combustion with chemical reactions.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
35 #include "psiReactionThermo.H"
36 #include "CombustionModel.H"
37 #include "multivariateScheme.H"
38 #include "pimpleControl.H"
39 #include "pressureControl.H"
40 #include "fvOptions.H"
41 #include "localEulerDdtScheme.H"
42 #include "fvcSmooth.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 int main(int argc, char *argv[])
47 {
48  #include "postProcess.H"
49 
50  #include "setRootCaseLists.H"
51  #include "createTime.H"
52  #include "createMesh.H"
53  #include "createControl.H"
54  #include "createTimeControls.H"
55  #include "initContinuityErrs.H"
56  #include "createFields.H"
57  #include "createFieldRefs.H"
58 
59  turbulence->validate();
60 
61  if (!LTS)
62  {
63  #include "compressibleCourantNo.H"
64  #include "setInitialDeltaT.H"
65  }
66 
67  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69  Info<< "\nStarting time loop\n" << endl;
70 
71  while (pimple.run(runTime))
72  {
73  #include "readTimeControls.H"
74 
75  if (LTS)
76  {
77  #include "setRDeltaT.H"
78  }
79  else
80  {
81  #include "compressibleCourantNo.H"
82  #include "setDeltaT.H"
83  }
84 
85  runTime++;
86 
87  Info<< "Time = " << runTime.timeName() << nl << endl;
88 
89  #include "rhoEqn.H"
90 
91  while (pimple.loop())
92  {
93  #include "UEqn.H"
94  #include "YEqn.H"
95  #include "EEqn.H"
96 
97  // --- Pressure corrector loop
98  while (pimple.correct())
99  {
100  #include "pEqn.H"
101  }
102 
103  if (pimple.turbCorr())
104  {
105  turbulence->correct();
106  thermophysicalTransport->correct();
107  }
108  }
109 
110  rho = thermo.rho();
111 
112  runTime.write();
113 
114  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
115  << " ClockTime = " << runTime.elapsedClockTime() << " s"
116  << nl << endl;
117  }
118 
119  Info<< "End\n" << endl;
120 
121  return 0;
122 }
123 
124 
125 // ************************************************************************* //
pimpleNoLoopControl & pimple
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Read the control parameters used by setDeltaT.
rhoReactionThermo & thermo
Definition: createFields.H:28
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
bool LTS
Definition: createRDeltaT.H:1
messageStream Info
Execute application functionObjects to post-process existing results.
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Read the control parameters used by setDeltaT.