multiphaseEulerFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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 momentum, heat and mass transfer.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "fvCFD.H"
37 #include "phaseSystem.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 "setRootCaseLists.H"
50  #include "createTime.H"
51  #include "createMesh.H"
52  #include "createControl.H"
53  #include "createTimeControls.H"
54  #include "createFields.H"
55  #include "createFieldRefs.H"
56 
57  if (!LTS)
58  {
59  #include "CourantNo.H"
60  #include "setInitialDeltaT.H"
61  }
62 
63  Switch faceMomentum
64  (
65  pimple.dict().lookupOrDefault<Switch>("faceMomentum", false)
66  );
67  Switch partialElimination
68  (
69  pimple.dict().lookupOrDefault<Switch>("partialElimination", false)
70  );
71 
72  #include "createRDeltaTf.H"
73 
74  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75 
76  Info<< "\nStarting time loop\n" << endl;
77 
78  while (pimple.run(runTime))
79  {
80  #include "readTimeControls.H"
81 
82  int nEnergyCorrectors
83  (
84  pimple.dict().lookupOrDefault<int>("nEnergyCorrectors", 1)
85  );
86 
87  if (LTS)
88  {
89  #include "setRDeltaT.H"
90  if (faceMomentum)
91  {
92  #include "setRDeltaTf.H"
93  }
94  }
95  else
96  {
97  #include "CourantNo.H"
98  #include "setDeltaT.H"
99  }
100 
101  runTime++;
102  Info<< "Time = " << runTime.timeName() << nl << endl;
103 
104  // --- Pressure-velocity PIMPLE corrector loop
105  while (pimple.loop())
106  {
107  fluid.solve(rAUs, rAUfs);
108  fluid.correct();
109  fluid.correctContinuityError();
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  }
119  else
120  {
121  #include "pU/UEqns.H"
122  #include "EEqns.H"
123  #include "pU/pEqn.H"
124  }
125 
126  fluid.correctKinematics();
127 
128  if (pimple.turbCorr())
129  {
130  fluid.correctTurbulence();
131  }
132  }
133 
134  runTime.write();
135 
136  Info<< "ExecutionTime = "
137  << runTime.elapsedCpuTime()
138  << " s\n\n" << endl;
139  }
140 
141  Info<< "End\n" << endl;
142 
143  return 0;
144 }
145 
146 
147 // ************************************************************************* //
pimpleNoLoopControl & pimple
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Read the control parameters used by setDeltaT.
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);label pRefCell=0;scalar pRefValue=0.0;if(fluid.incompressible()){ p=max(p_rgh+fluid.rho() *gh, pMin);if(p_rgh.needReference()) { setRefCell(p, p_rgh, pimple.dict(), pRefCell, pRefValue);p+=dimensionedScalar("p", p.dimensions(), pRefValue - getRefCellValue(p, pRefCell));p_rgh=p - fluid.rho() *gh;}}mesh.setFluxRequired(p_rgh.name());PtrList< volScalarField > rAUs
Definition: createFields.H:67
phaseSystem & fluid
Definition: createFields.H:11
static const char nl
Definition: Ostream.H:260
bool LTS
Definition: createRDeltaT.H:1
PtrList< surfaceScalarField > rAUfs
Definition: createFields.H:68
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.