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