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-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  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"
39 #include "pimpleControl.H"
40 #include "localEulerDdtScheme.H"
41 #include "fvcSmooth.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 int main(int argc, char *argv[])
46 {
47  #include "postProcess.H"
48 
49  #include "setRootCase.H"
50  #include "createTime.H"
51  #include "createMesh.H"
52  #include "createControl.H"
53  #include "createTimeControls.H"
54  #include "createRDeltaT.H"
55  #include "createFields.H"
56  #include "createFieldRefs.H"
57 
58  if (!LTS)
59  {
60  #include "CourantNo.H"
61  #include "setInitialDeltaT.H"
62  }
63 
64  Switch faceMomentum
65  (
66  pimple.dict().lookupOrDefault<Switch>("faceMomentum", false)
67  );
68 
69  Switch implicitPhasePressure
70  (
71  mesh.solverDict(alpha1.name()).lookupOrDefault<Switch>
72  (
73  "implicitPhasePressure", false
74  )
75  );
76 
77  #include "pUf/createDDtU.H"
78 
79  int nEnergyCorrectors
80  (
81  pimple.dict().lookupOrDefault<int>("nEnergyCorrectors", 1)
82  );
83 
84  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85 
86  Info<< "\nStarting time loop\n" << endl;
87 
88  while (runTime.run())
89  {
90  #include "readTimeControls.H"
91 
92  if (LTS)
93  {
94  #include "setRDeltaT.H"
95  }
96  else
97  {
98  #include "CourantNos.H"
99  #include "setDeltaT.H"
100  }
101 
102  runTime++;
103  Info<< "Time = " << runTime.timeName() << nl << endl;
104 
105  // --- Pressure-velocity PIMPLE corrector loop
106  while (pimple.loop())
107  {
108  fluid.solve();
109  fluid.correct();
110 
111  #include "YEqns.H"
112 
113  if (faceMomentum)
114  {
115  #include "pUf/UEqns.H"
116  #include "EEqns.H"
117  #include "pUf/pEqn.H"
118  #include "pUf/DDtU.H"
119  }
120  else
121  {
122  #include "pU/UEqns.H"
123  #include "EEqns.H"
124  #include "pU/pEqn.H"
125  }
126 
127  fluid.correctKinematics();
128 
129  if (pimple.turbCorr())
130  {
131  fluid.correctTurbulence();
132  }
133  }
134 
135  runTime.write();
136 
137  Info<< "ExecutionTime = "
138  << runTime.elapsedCpuTime()
139  << " s\n\n" << endl;
140  }
141 
142  Info<< "End\n" << endl;
143 
144  return 0;
145 }
146 
147 
148 // ************************************************************************* //
multiphaseSystem & fluid
Definition: createFields.H:10
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Read the control parameters used by setDeltaT.
dynamicFvMesh & mesh
static const char nl
Definition: Ostream.H:262
bool LTS
Definition: createRDeltaT.H:1
volScalarField & alpha1
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.