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 | 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  interMixingFoam
26 
27 Description
28  Solver for 3 incompressible fluids, two of which are miscible,
29  using a VOF method to capture the interface.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
34 #include "CMULES.H"
35 #include "subCycle.H"
38 #include "pimpleControl.H"
39 #include "fvOptions.H"
40 #include "CorrectPhi.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 "setRootCase.H"
51  #include "createTime.H"
52  #include "createMesh.H"
53  #include "createControl.H"
54  #include "createTimeControls.H"
55  #include "createRDeltaT.H"
56  #include "initContinuityErrs.H"
57  #include "createFields.H"
58  #include "createFvOptions.H"
59  #include "correctPhi.H"
60 
61  turbulence->validate();
62 
63  if (!LTS)
64  {
65  #include "readTimeControls.H"
66  #include "CourantNo.H"
67  #include "setInitialDeltaT.H"
68  }
69 
70  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72  Info<< "\nStarting time loop\n" << endl;
73 
74  while (runTime.run())
75  {
76  #include "readTimeControls.H"
77 
78  if (LTS)
79  {
80  #include "setRDeltaT.H"
81  }
82  else
83  {
84  #include "CourantNo.H"
85  #include "alphaCourantNo.H"
86  #include "setDeltaT.H"
87  }
88 
89  runTime++;
90 
91  Info<< "Time = " << runTime.timeName() << nl << endl;
92 
93  // --- Pressure-velocity PIMPLE corrector loop
94  while (pimple.loop())
95  {
96  #include "alphaControls.H"
97  #include "alphaEqnsSubCycle.H"
98 
99  mixture.correct();
100 
101  #include "UEqn.H"
102 
103  // --- Pressure corrector loop
104  while (pimple.correct())
105  {
106  #include "pEqn.H"
107  }
108 
109  if (pimple.turbCorr())
110  {
111  turbulence->correct();
112  }
113  }
114 
115  #include "continuityErrs.H"
116 
117  runTime.write();
118 
119  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
120  << " ClockTime = " << runTime.elapsedClockTime() << " s"
121  << nl << endl;
122  }
123 
124  Info<< "End\n" << endl;
125 
126  return 0;
127 }
128 
129 
130 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Read the control parameters used by setDeltaT.
static const char nl
Definition: Ostream.H:262
bool LTS
Definition: createRDeltaT.H:1
Info<< "Reading field p_rgh\n"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating phaseChangeTwoPhaseMixture\n"<< endl;autoPtr< phaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:33
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.
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.