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