reactingMultiphaseEulerFoam.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  reactingMultiphaseEulerFoam
26 
27 Description
28  Solver for a system of any number of compressible fluid phases with a
29  common pressure, but otherwise separate properties. The type of phase model
30  is run time selectable and can optionally represent multiple species and
31  in-phase reactions. The phase system is also run time selectable and can
32  optionally represent different types of momentun, heat and mass transfer.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "fvCFD.H"
37 #include "multiphaseSystem.H"
38 #include "pimpleControl.H"
39 #include "localEulerDdtScheme.H"
40 #include "fvcSmooth.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 int main(int argc, char *argv[])
45 {
46  #include "postProcess.H"
47 
48  #include "setRootCase.H"
49  #include "createTime.H"
50  #include "createMesh.H"
51  #include "createControl.H"
52  #include "createTimeControls.H"
53  #include "createFields.H"
54  #include "createFieldRefs.H"
55 
56  if (!LTS)
57  {
58  #include "CourantNo.H"
59  #include "setInitialDeltaT.H"
60  }
61 
62  // Switch faceMomentum
63  // (
64  // pimple.dict().lookupOrDefault<Switch>("faceMomentum", false)
65  // );
66 
67  // Switch implicitPhasePressure
68  // (
69  // mesh.solverDict(alpha1.name()).lookupOrDefault<Switch>
70  // (
71  // "implicitPhasePressure", false
72  // )
73  // );
74 
75  //#include "pUf/createDDtU.H"
76 
77  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78 
79  Info<< "\nStarting time loop\n" << endl;
80 
81  while (runTime.run())
82  {
83  #include "readTimeControls.H"
84 
85  int nEnergyCorrectors
86  (
87  pimple.dict().lookupOrDefault<int>("nEnergyCorrectors", 1)
88  );
89 
90  if (LTS)
91  {
92  #include "setRDeltaT.H"
93  }
94  else
95  {
96  #include "CourantNo.H"
97  #include "setDeltaT.H"
98  }
99 
100  runTime++;
101  Info<< "Time = " << runTime.timeName() << nl << endl;
102 
103  // --- Pressure-velocity PIMPLE corrector loop
104  while (pimple.loop())
105  {
106  fluid.solve();
107  fluid.correct();
108 
109  #include "YEqns.H"
110 
111  // if (faceMomentum)
112  // {
113  // #include "pUf/UEqns.H"
114  // #include "EEqns.H"
115  // #include "pUf/pEqn.H"
116  // #include "pUf/DDtU.H"
117  // }
118  // else
119  {
120  #include "pU/UEqns.H"
121  #include "EEqns.H"
122  #include "pU/pEqn.H"
123  }
124 
125  fluid.correctKinematics();
126 
127  if (pimple.turbCorr())
128  {
129  fluid.correctTurbulence();
130  }
131  }
132 
133  runTime.write();
134 
135  Info<< "ExecutionTime = "
136  << runTime.elapsedCpuTime()
137  << " s\n\n" << endl;
138  }
139 
140  Info<< "End\n" << endl;
141 
142  return 0;
143 }
144 
145 
146 // ************************************************************************* //
multiphaseSystem & fluid
Definition: createFields.H:11
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Read the control parameters used by setDeltaT.
static const char nl
Definition: Ostream.H:262
bool LTS
Definition: createRDeltaT.H:1
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.