reactingParcelFoam.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-2015 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  reactingParcelFoam
26 
27 Description
28  Transient PIMPLE solver for compressible, laminar or turbulent flow with
29  reacting multiphase Lagrangian parcels, including run-time selectable
30  finite volume options, e.g. sources, constraints
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #include "fvCFD.H"
37 #include "rhoCombustionModel.H"
38 #include "radiationModel.H"
39 #include "fvIOoptionList.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 "setRootCase.H"
50 
51  #include "createTime.H"
52  #include "createMesh.H"
53  #include "readGravitationalAcceleration.H"
54 
55  pimpleControl pimple(mesh);
56 
57  #include "createTimeControls.H"
58  #include "createRDeltaT.H"
59  #include "initContinuityErrs.H"
60  #include "createFields.H"
61  #include "createRadiationModel.H"
62  #include "createClouds.H"
63  #include "createMRF.H"
64  #include "createFvOptions.H"
65 
66  if (!LTS)
67  {
68  #include "compressibleCourantNo.H"
69  #include "setInitialDeltaT.H"
70  }
71 
72  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
74  Info<< "\nStarting time loop\n" << endl;
75 
76  while (runTime.run())
77  {
78  #include "readTimeControls.H"
79 
80  if (!LTS)
81  {
82  #include "compressibleCourantNo.H"
83  #include "setDeltaT.H"
84  }
85 
86  runTime++;
87 
88  Info<< "Time = " << runTime.timeName() << nl << endl;
89 
90  parcels.evolve();
91 
92  if (LTS)
93  {
94  #include "setRDeltaT.H"
95  }
96 
97  #include "rhoEqn.H"
98 
99  // --- Pressure-velocity PIMPLE corrector loop
100  while (pimple.loop())
101  {
102  #include "UEqn.H"
103  #include "YEqn.H"
104  #include "EEqn.H"
105 
106  // --- Pressure corrector loop
107  while (pimple.correct())
108  {
109  #include "pEqn.H"
110  }
111 
112  if (pimple.turbCorr())
113  {
114  turbulence->correct();
115  }
116  }
117 
118  rho = thermo.rho();
119 
120  runTime.write();
121 
122  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
123  << " ClockTime = " << runTime.elapsedClockTime() << " s"
124  << nl << endl;
125  }
126 
127  Info<< "End\n" << endl;
128 
129  return 0;
130 }
131 
132 
133 // ************************************************************************* //
bool LTS
Definition: createRDeltaT.H:1
messageStream Info
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
dynamicFvMesh & mesh
const dictionary & pimple
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Read the control parameters used by setDeltaT.
int main(int argc, char *argv[])
Definition: postCalc.C:54
psiReactionThermo & thermo
Definition: createFields.H:32
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Read the control parameters used by setDeltaT.