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 | 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  coalChemistryFoam
28 
29 Group
30  grpLagrangianSolvers
31 
32 Description
33  Transient solver for compressible, turbulent flow, with coal and limestone
34  particle clouds, an energy source, and combustion.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "fvCFD.H"
40 #include "basicThermoCloud.H"
41 #include "coalCloud.H"
42 #include "psiReactionThermo.H"
43 #include "CombustionModel.H"
44 #include "fvOptions.H"
45 #include "radiationModel.H"
46 #include "SLGThermo.H"
47 #include "pimpleControl.H"
48 #include "pressureControl.H"
49 #include "localEulerDdtScheme.H"
50 #include "fvcSmooth.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 int main(int argc, char *argv[])
55 {
56  argList::addNote
57  (
58  "Transient solver for compressible, turbulent flow"
59  " with coal and limestone clouds, energy sources and combustion."
60  );
61 
62  #include "postProcess.H"
63 
64  #include "addCheckCaseOptions.H"
65  #include "setRootCaseLists.H"
66  #include "createTime.H"
67  #include "createMesh.H"
68  #include "createControl.H"
69  #include "createTimeControls.H"
70  #include "createFields.H"
71  #include "createFieldRefs.H"
72  #include "initContinuityErrs.H"
73 
74  turbulence->validate();
75 
76  if (!LTS)
77  {
78  #include "compressibleCourantNo.H"
79  #include "setInitialDeltaT.H"
80  }
81 
82  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83 
84  Info<< "\nStarting time loop\n" << endl;
85 
86  while (runTime.run())
87  {
88  #include "readTimeControls.H"
89 
90  if (LTS)
91  {
92  #include "setRDeltaT.H"
93  }
94  else
95  {
96  #include "compressibleCourantNo.H"
97  #include "setDeltaT.H"
98  }
99 
100  ++runTime;
101 
102  Info<< "Time = " << runTime.timeName() << nl << endl;
103 
104  rhoEffLagrangian = coalParcels.rhoEff() + limestoneParcels.rhoEff();
105  pDyn = 0.5*rho*magSqr(U);
106 
107  coalParcels.evolve();
108 
109  limestoneParcels.evolve();
110 
111  #include "rhoEqn.H"
112 
113  // --- Pressure-velocity PIMPLE corrector loop
114  while (pimple.loop())
115  {
116  #include "UEqn.H"
117  #include "YEqn.H"
118  #include "EEqn.H"
119 
120  // --- Pressure corrector loop
121  while (pimple.correct())
122  {
123  #include "pEqn.H"
124  }
125 
126  if (pimple.turbCorr())
127  {
128  turbulence->correct();
129  }
130  }
131 
132  rho = thermo.rho();
133 
134  runTime.write();
135 
136  runTime.printExecutionTime(Info);
137  }
138 
139  Info<< "End\n" << endl;
140 
141  return 0;
142 }
143 
144 
145 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
fvOptions.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
pDyn
volScalarField pDyn(IOobject("pDyn", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimPressure, Zero))
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
SLGThermo.H
rho
rho
Definition: readInitialConditions.H:88
pimpleControl.H
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
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)
CombustionModel.H
rhoEffLagrangian
volScalarField rhoEffLagrangian(IOobject("rhoEffLagrangian", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimDensity, Zero))
coalCloud.H
LTS
bool LTS
Definition: createRDeltaT.H:1
postProcess.H
Execute application functionObjects to post-process existing results.
pressureControl.H
U
U
Definition: pEqn.H:72
basicThermoCloud.H
psiReactionThermo.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
radiationModel.H
turbulentFluidThermoModel.H
fvcSmooth.H
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...