interMixingFoam.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) 2011-2017 OpenFOAM Foundation
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  interMixingFoam
28 
29 Group
30  grpMultiphaseSolvers
31 
32 Description
33  Solver for 3 incompressible fluids, two of which are miscible, using a VOF
34  method to capture the interface, with optional mesh motion and mesh topology
35  changes including adaptive re-meshing.
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "fvCFD.H"
40 #include "dynamicFvMesh.H"
41 #include "CMULES.H"
42 #include "localEulerDdtScheme.H"
43 #include "subCycle.H"
46 #include "pimpleControl.H"
47 #include "fvOptions.H"
48 #include "CorrectPhi.H"
49 #include "fvcSmooth.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 int main(int argc, char *argv[])
54 {
55  argList::addNote
56  (
57  "Solver for three incompressible fluids (two of which are immiscible)"
58  " using VOF phase-fraction based interface capturing.\n"
59  "With optional mesh motion and mesh topology changes including"
60  " adaptive re-meshing."
61  );
62 
63  #include "postProcess.H"
64 
65  #include "addCheckCaseOptions.H"
66  #include "setRootCaseLists.H"
67  #include "createTime.H"
68  #include "createDynamicFvMesh.H"
69  #include "initContinuityErrs.H"
70  #include "createDyMControls.H"
71  #include "createFields.H"
72  #include "initCorrectPhi.H"
73  #include "createUfIfPresent.H"
74 
75  turbulence->validate();
76 
77  if (!LTS)
78  {
79  #include "readTimeControls.H"
80  #include "CourantNo.H"
81  #include "setInitialDeltaT.H"
82  }
83 
84  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85 
86  Info<< "\nStarting time loop\n" << endl;
87 
88  while (runTime.run())
89  {
90  #include "readDyMControls.H"
91 
92  if (LTS)
93  {
94  #include "setRDeltaT.H"
95  }
96  else
97  {
98  #include "CourantNo.H"
99  #include "alphaCourantNo.H"
100  #include "setDeltaT.H"
101  }
102 
103  ++runTime;
104 
105  Info<< "Time = " << runTime.timeName() << nl << endl;
106 
107  // --- Pressure-velocity PIMPLE corrector loop
108  while (pimple.loop())
109  {
110  if (pimple.firstIter() || moveMeshOuterCorrectors)
111  {
112  mesh.update();
113 
114  if (mesh.changing())
115  {
116  gh = (g & mesh.C()) - ghRef;
117  ghf = (g & mesh.Cf()) - ghRef;
118 
119  MRF.update();
120 
121  if (correctPhi)
122  {
123  // Calculate absolute flux
124  // from the mapped surface velocity
125  phi = mesh.Sf() & Uf();
126 
127  #include "correctPhi.H"
128 
129  // Make the flux relative to the mesh motion
131 
132  mixture.correct();
133  }
134 
135  if (checkMeshCourantNo)
136  {
137  #include "meshCourantNo.H"
138  }
139  }
140  }
141 
142  #include "alphaControls.H"
143  #include "alphaEqnSubCycle.H"
144 
145  mixture.correct();
146 
147  #include "UEqn.H"
148 
149  // --- Pressure corrector loop
150  while (pimple.correct())
151  {
152  #include "pEqn.H"
153  }
154 
155  if (pimple.turbCorr())
156  {
157  turbulence->correct();
158  }
159  }
160 
161  #include "continuityErrs.H"
162 
163  runTime.write();
164 
165  runTime.printExecutionTime(Info);
166  }
167 
168  Info<< "End\n" << endl;
169 
170  return 0;
171 }
172 
173 
174 // ************************************************************************* //
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Read the control parameters used by setDeltaT.
correctPhi
checkMeshCourantNo
const surfaceScalarField & ghf
dynamicFvMesh & mesh
IOMRFZoneList & MRF
autoPtr< surfaceVectorField > Uf
const uniformDimensionedVectorField & g
bool LTS
Definition: createRDeltaT.H:1
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
Calculates and outputs the mean and maximum Courant Numbers.
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
Required Classes.
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Creates and initialises the velocity field Uf if required.