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-2022 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  multiphaseEulerFoam
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 "pressureReference.H"
41 #include "localEulerDdtScheme.H"
42 #include "fvcSmooth.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 int main(int argc, char *argv[])
47 {
48  #include "postProcess.H"
49 
50  #include "setRootCaseLists.H"
51  #include "createTime.H"
52  #include "createMesh.H"
53  #include "createDyMControls.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 "readDyMControls.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 
102 
103  // Store divU from the previous mesh so that it can be
104  // mapped and used in correctPhi to ensure the corrected phi
105  // has the same divergence
106  tmp<volScalarField> divU;
107 
108  if (correctPhi && mesh.topoChanged())
109  {
110  // Construct and register divU for mapping
111  divU = new volScalarField
112  (
113  "divU0",
114  fvc::div
115  (
116  fvc::absolute(phi, fluid.movingPhases()[0].U())
117  )
118  );
119  }
120 
121  mesh.update();
122 
123  runTime++;
124 
125  Info<< "Time = " << runTime.userTimeName() << nl << endl;
126 
127  // --- Pressure-velocity PIMPLE corrector loop
128  while (pimple.loop())
129  {
130  if (!pimple.flow())
131  {
132  if (pimple.models())
133  {
134  fvModels.correct();
135  }
136 
137  if (pimple.thermophysics())
138  {
139  fluid.solve(rAUs, rAUfs);
140  fluid.correct();
141  fluid.correctContinuityError();
142 
143  #include "YEqns.H"
144  #include "EEqns.H"
145  #include "pEqnComps.H"
146 
148  {
149  phases[phasei].divU(-pEqnComps[phasei] & p_rgh);
150  }
151  }
152  }
153  else
154  {
155  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
156  {
157  if (correctPhi && !divU.valid())
158  {
159  // Construct and register divU for mapping
160  divU = new volScalarField
161  (
162  "divU0",
163  fvc::div
164  (
165  fvc::absolute(phi, fluid.movingPhases()[0].U())
166  )
167  );
168  }
169 
170  // Move the mesh
171  mesh.move();
172 
173  if (mesh.changing())
174  {
175  gh = (g & mesh.C()) - ghRef;
176  ghf = (g & mesh.Cf()) - ghRef;
177 
178  fluid.meshUpdate();
179 
180  if (correctPhi)
181  {
182  fluid.correctPhi
183  (
184  p_rgh,
185  divU,
187  pimple
188  );
189  }
190 
191  if (checkMeshCourantNo)
192  {
193  #include "meshCourantNo.H"
194  }
195  }
196 
197  divU.clear();
198  }
199 
200  if (pimple.models())
201  {
202  fvModels.correct();
203  }
204 
205  fluid.solve(rAUs, rAUfs);
206  fluid.correct();
207  fluid.correctContinuityError();
208 
209  if (pimple.thermophysics())
210  {
211  #include "YEqns.H"
212  }
213 
214  if (faceMomentum)
215  {
216  #include "pUf/UEqns.H"
217 
218  if (pimple.thermophysics())
219  {
220  #include "EEqns.H"
221  }
222 
223  #include "pUf/pEqn.H"
224  }
225  else
226  {
227  #include "pU/UEqns.H"
228 
229  if (pimple.thermophysics())
230  {
231  #include "EEqns.H"
232  }
233 
234  #include "pU/pEqn.H"
235  }
236 
237  fluid.correctKinematics();
238 
239  if (pimple.turbCorr())
240  {
241  fluid.correctTurbulence();
242  }
243  }
244  }
245 
246  runTime.write();
247 
248  Info<< "ExecutionTime = "
249  << runTime.elapsedCpuTime()
250  << " s\n\n" << endl;
251  }
252 
253  Info<< "End\n" << endl;
254 
255  return 0;
256 }
257 
258 
259 // ************************************************************************* //
pressureReference & pressureReference
volScalarField & p_rgh
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
pimpleNoLoopControl & pimple
label phasei
Definition: pEqn.H:27
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:353
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:255
correctPhi
checkMeshCourantNo
fvMesh & mesh
PtrList< fvScalarMatrix > pEqnComps(phases.size())
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);if(fluid.found("pMin")){ IOWarningInFunction(fluid)<< "Pressure limits, pMin and pMax, are now read from "<< pimple.dict().name()<< endl;}pressureReference pressureReference(p, p_rgh, pimple.dict(), fluid.incompressible());if(fluid.incompressible()){ p=p_rgh+fluid.rho() *gh;}if(p_rgh.needReference() &&fluid.incompressible()){ p+=dimensionedScalar("p", p.dimensions(), pressureReference.refValue() - getRefCellValue(p, pressureReference.refCell()));}p_rgh=p - fluid.rho() *gh;mesh.schemes().setFluxRequired(p_rgh.name());PtrList< volScalarField > rAUs
Definition: createFields.H:67
phaseSystem & fluid
Definition: createFields.H:11
phi
Definition: correctPhi.H:3
const surfaceScalarField & ghf
static const char nl
Definition: Ostream.H:260
bool LTS
Definition: createRDeltaT.H:1
PtrList< surfaceScalarField > rAUfs
Definition: createFields.H:68
moveMeshOuterCorrectors
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
Foam::fvModels & fvModels
Calculates and outputs the mean and maximum Courant Numbers.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
phaseSystem::phaseModelList & phases
Definition: createFields.H:12
const volScalarField & gh
messageStream Info
const dimensionedVector & g
Execute application functionObjects to post-process existing results.
Provides functions smooth spread and sweep which use the FvFaceCellWave algorithm to smooth and redis...