coalChemistryFoam.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 | Copyright (C) 2011-2016 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  coalChemistryFoam
26 
27 Description
28  Transient solver for compressible, turbulent flow, with coal and limestone
29  particle clouds, an energy source, and combustion.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
35 #include "basicThermoCloud.H"
36 #include "coalCloud.H"
37 #include "psiCombustionModel.H"
38 #include "fvOptions.H"
39 #include "radiationModel.H"
40 #include "SLGThermo.H"
41 #include "pimpleControl.H"
42 #include "localEulerDdtScheme.H"
43 #include "fvcSmooth.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 int main(int argc, char *argv[])
48 {
49  #include "postProcess.H"
50 
51  #include "setRootCase.H"
52  #include "createTime.H"
53  #include "createMesh.H"
54  #include "createControl.H"
55  #include "createTimeControls.H"
56  #include "createRDeltaT.H"
57  #include "createFields.H"
58  #include "createFvOptions.H"
59  #include "initContinuityErrs.H"
60 
61  turbulence->validate();
62 
63  if (!LTS)
64  {
65  #include "compressibleCourantNo.H"
66  #include "setInitialDeltaT.H"
67  }
68 
69  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71  Info<< "\nStarting time loop\n" << endl;
72 
73  while (runTime.run())
74  {
75  #include "readTimeControls.H"
76 
77  if (LTS)
78  {
79  #include "setRDeltaT.H"
80  }
81  else
82  {
83  #include "compressibleCourantNo.H"
84  #include "setDeltaT.H"
85  }
86 
87  runTime++;
88 
89  Info<< "Time = " << runTime.timeName() << nl << endl;
90 
91  rhoEffLagrangian = coalParcels.rhoEff() + limestoneParcels.rhoEff();
92  pDyn = 0.5*rho*magSqr(U);
93 
94  coalParcels.evolve();
95 
96  limestoneParcels.evolve();
97 
98  #include "rhoEqn.H"
99 
100  // --- Pressure-velocity PIMPLE corrector loop
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  #include "pEqn.H"
111  }
112 
113  if (pimple.turbCorr())
114  {
115  turbulence->correct();
116  }
117  }
118 
119  rho = thermo.rho();
120 
121  runTime.write();
122 
123  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
124  << " ClockTime = " << runTime.elapsedClockTime() << " s"
125  << nl << endl;
126  }
127 
128  Info<< "End\n" << endl;
129 
130  return 0;
131 }
132 
133 
134 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
U
Definition: pEqn.H:83
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Read the control parameters used by setDeltaT.
psiReactionThermo & thermo
Definition: createFields.H:31
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:262
bool LTS
Definition: createRDeltaT.H:1
volScalarField rhoEffLagrangian(IOobject("rhoEffLagrangian", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar("zero", dimDensity, 0.0))
messageStream Info
volScalarField pDyn(IOobject("pDyn", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar("zero", dimPressure, 0.0))
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.