icoReactingMultiphaseInterFoam.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2017-2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  icoReactingMultiphaseInterFoam
28 
29 Group
30  grpMultiphaseSolvers
31 
32 Description
33  Solver for N incompressible, non-isothermal immiscible fluids with
34  phase-change. Uses a VOF (volume of fluid) phase-fraction based interface
35  capturing approach.
36 
37  The momentum, energy and other fluid properties are of the "mixture" and a
38  single momentum equation is solved.
39 
40  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #include "fvCFD.H"
45 #include "dynamicFvMesh.H"
46 #include "subCycle.H"
47 #include "multiphaseSystem.H"
49 #include "pimpleControl.H"
50 #include "fvOptions.H"
51 #include "radiationModel.H"
52 #include "CorrectPhi.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 int main(int argc, char *argv[])
57 {
58  argList::addNote
59  (
60  "Solver for N incompressible, non-isothermal immiscible fluids with"
61  " phase-change,"
62  " using VOF phase-fraction based interface capturing.\n"
63  "With optional mesh motion and mesh topology changes including"
64  " adaptive re-meshing."
65  );
66 
67  #include "postProcess.H"
68 
69  #include "setRootCaseLists.H"
70  #include "createTime.H"
71  #include "createDynamicFvMesh.H"
72 
73  #include "initContinuityErrs.H"
74  #include "createDyMControls.H"
75 
76  #include "createFields.H"
77  #include "createFieldRefs.H"
78 
79  #include "initCorrectPhi.H"
80  #include "createUfIfPresent.H"
81 
82  #include "createFvOptions.H"
83  #include "CourantNo.H"
84  #include "setInitialDeltaT.H"
85 
86  turbulence->validate();
87 
88  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
89 
90  Info<< "\nStarting time loop\n" << endl;
91 
92  while (runTime.run())
93  {
94  #include "readDyMControls.H"
95 
96  #include "CourantNo.H"
97  #include "alphaCourantNo.H"
98  #include "setDeltaT.H"
99 
100  ++runTime;
101 
102  Info<< "Time = " << runTime.timeName() << nl << endl;
103 
104  fluid.correctMassSources(T);
105  fluid.solve();
106  rho = fluid.rho();
107 
108  // --- Pressure-velocity PIMPLE corrector loop
109  while (pimple.loop())
110  {
111  if (pimple.firstIter() || moveMeshOuterCorrectors)
112  {
113  mesh.update();
114 
115  if (mesh.changing())
116  {
117  gh = (g & mesh.C()) - ghRef;
118  ghf = (g & mesh.Cf()) - ghRef;
119 
120  if (correctPhi)
121  {
122  // Calculate absolute flux
123  // from the mapped surface velocity
124  phi = mesh.Sf() & Uf();
125 
126  #include "correctPhi.H"
127 
128  // Make the flux relative to the mesh motion
130  }
131  }
132  }
133 
134  #include "YEqns.H"
135  #include "TEqn.H"
136 
137  if (pimple.frozenFlow())
138  {
139  continue;
140  }
141 
142  #include "UEqn.H"
143 
144  // --- Pressure corrector loop
145  while (pimple.correct())
146  {
147  #include "pEqn.H"
148  }
149 
150  if (pimple.turbCorr())
151  {
152  turbulence->correct();
153  }
154  }
155 
156  rho = fluid.rho();
157 
158  runTime.write();
159 
160  runTime.printExecutionTime(Info);
161  }
162 
163  Info<< "End\n" << endl;
164 
165  return 0;
166 }
167 
168 
169 // ************************************************************************* //
twoPhaseSystem & fluid
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
correctPhi
const surfaceScalarField & ghf
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
const uniformDimensionedVectorField & g
pimpleControl & pimple
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
moveMeshOuterCorrectors
const volScalarField & T
U
Definition: pEqn.H:72
const volScalarField & gh
messageStream Info
Information stream (stdout output on master, null elsewhere)
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:70
Creates and initialises the velocity field Uf if required.