interFoam.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  interFoam
26 
27 Description
28  Solver for 2 incompressible, isothermal immiscible fluids using a VOF
29  (volume of fluid) phase-fraction based interface capturing approach.
30 
31  The momentum and other fluid properties are of the "mixture" and a single
32  momentum equation is solved.
33 
34  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
35 
36  For a two-fluid approach see twoPhaseEulerFoam.
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #include "fvCFD.H"
41 #include "CMULES.H"
42 #include "EulerDdtScheme.H"
43 #include "localEulerDdtScheme.H"
44 #include "CrankNicolsonDdtScheme.H"
45 #include "subCycle.H"
48 #include "pimpleControl.H"
49 #include "fvOptions.H"
50 #include "CorrectPhi.H"
51 #include "localEulerDdtScheme.H"
52 #include "fvcSmooth.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 int main(int argc, char *argv[])
57 {
58  #include "postProcess.H"
59 
60  #include "setRootCase.H"
61  #include "createTime.H"
62  #include "createMesh.H"
63  #include "createControl.H"
64  #include "createTimeControls.H"
65  #include "createRDeltaT.H"
66  #include "initContinuityErrs.H"
67  #include "createFields.H"
68  #include "createFvOptions.H"
69  #include "correctPhi.H"
70 
71  turbulence->validate();
72 
73  if (!LTS)
74  {
75  #include "readTimeControls.H"
76  #include "CourantNo.H"
77  #include "setInitialDeltaT.H"
78  }
79 
80  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81 
82  Info<< "\nStarting time loop\n" << endl;
83 
84  while (runTime.run())
85  {
86  #include "readTimeControls.H"
87 
88  if (LTS)
89  {
90  #include "setRDeltaT.H"
91  }
92  else
93  {
94  #include "CourantNo.H"
95  #include "alphaCourantNo.H"
96  #include "setDeltaT.H"
97  }
98 
99  runTime++;
100 
101  Info<< "Time = " << runTime.timeName() << nl << endl;
102 
103  // --- Pressure-velocity PIMPLE corrector loop
104  while (pimple.loop())
105  {
106  #include "alphaControls.H"
107  #include "alphaEqnSubCycle.H"
108 
109  mixture.correct();
110 
111  #include "UEqn.H"
112 
113  // --- Pressure corrector loop
114  while (pimple.correct())
115  {
116  #include "pEqn.H"
117  }
118 
119  if (pimple.turbCorr())
120  {
121  turbulence->correct();
122  }
123  }
124 
125  runTime.write();
126 
127  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
128  << " ClockTime = " << runTime.elapsedClockTime() << " s"
129  << nl << endl;
130  }
131 
132  Info<< "End\n" << endl;
133 
134  return 0;
135 }
136 
137 
138 // ************************************************************************* //
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.