coalChemistryFoam.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-2018 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 "psiReactionThermo.H"
38 #include "CombustionModel.H"
39 #include "fvOptions.H"
40 #include "radiationModel.H"
41 #include "SLGThermo.H"
42 #include "pimpleControl.H"
43 #include "pressureControl.H"
44 #include "localEulerDdtScheme.H"
45 #include "fvcSmooth.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 int main(int argc, char *argv[])
50 {
51  #include "postProcess.H"
52 
53  #include "setRootCaseLists.H"
54  #include "createTime.H"
55  #include "createMesh.H"
56  #include "createControl.H"
57  #include "createTimeControls.H"
58  #include "createFields.H"
59  #include "createFieldRefs.H"
60  #include "initContinuityErrs.H"
61 
62  turbulence->validate();
63 
64  if (!LTS)
65  {
66  #include "compressibleCourantNo.H"
67  #include "setInitialDeltaT.H"
68  }
69 
70  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72  Info<< "\nStarting time loop\n" << endl;
73 
74  while (runTime.run())
75  {
76  #include "readTimeControls.H"
77 
78  if (LTS)
79  {
80  #include "setRDeltaT.H"
81  }
82  else
83  {
84  #include "compressibleCourantNo.H"
85  #include "setDeltaT.H"
86  }
87 
88  runTime++;
89 
90  Info<< "Time = " << runTime.timeName() << nl << endl;
91 
92  rhoEffLagrangian = coalParcels.rhoEff() + limestoneParcels.rhoEff();
93  pDyn = 0.5*rho*magSqr(U);
94 
95  coalParcels.evolve();
96 
97  limestoneParcels.evolve();
98 
99  #include "rhoEqn.H"
100 
101  // --- Pressure-velocity PIMPLE corrector loop
102  while (pimple.loop())
103  {
104  #include "UEqn.H"
105  #include "YEqn.H"
106  #include "EEqn.H"
107 
108  // --- Pressure corrector loop
109  while (pimple.correct())
110  {
111  #include "pEqn.H"
112  }
113 
114  if (pimple.turbCorr())
115  {
116  turbulence->correct();
117  }
118  }
119 
120  rho = thermo.rho();
121 
122  runTime.write();
123 
124  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
125  << " ClockTime = " << runTime.elapsedClockTime() << " s"
126  << nl << endl;
127  }
128 
129  Info<< "End\n" << endl;
130 
131  return 0;
132 }
133 
134 
135 // ************************************************************************* //
pimpleNoLoopControl & pimple
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Read the control parameters used by setDeltaT.
rhoReactionThermo & thermo
Definition: createFields.H:28
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:265
bool LTS
Definition: createRDeltaT.H:1
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("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\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
volScalarField rhoEffLagrangian(IOobject("rhoEffLagrangian", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar("zero", dimDensity, 0.0))
U
Definition: pEqn.H:72
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.