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