compressibleInterFoam.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  compressibleInterFoam
26 
27 Description
28  Solver for 2 compressible, non-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 \*---------------------------------------------------------------------------*/
37 
38 #include "fvCFD.H"
39 #include "MULES.H"
40 #include "subCycle.H"
41 #include "rhoThermo.H"
42 #include "interfaceProperties.H"
43 #include "twoPhaseMixture.H"
44 #include "twoPhaseMixtureThermo.H"
46 #include "pimpleControl.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 int main(int argc, char *argv[])
51 {
52  #include "postProcess.H"
53 
54  #include "setRootCase.H"
55  #include "createTime.H"
56  #include "createMesh.H"
57  #include "createControl.H"
58  #include "createTimeControls.H"
59  #include "createFields.H"
60  #include "CourantNo.H"
61  #include "setInitialDeltaT.H"
62 
63  volScalarField& p = mixture.p();
64  volScalarField& T = mixture.T();
65  const volScalarField& psi1 = mixture.thermo1().psi();
66  const volScalarField& psi2 = mixture.thermo2().psi();
67 
68  turbulence->validate();
69 
70  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72  Info<< "\nStarting time loop\n" << endl;
73 
74  while (runTime.run())
75  {
76  #include "readTimeControls.H"
77  #include "CourantNo.H"
78  #include "setDeltaT.H"
79 
80  runTime++;
81 
82  Info<< "Time = " << runTime.timeName() << nl << endl;
83 
84  // --- Pressure-velocity PIMPLE corrector loop
85  while (pimple.loop())
86  {
87  #include "alphaEqnsSubCycle.H"
88 
89  // correct interface on first PIMPLE corrector
90  if (pimple.corr() == 1)
91  {
92  interface.correct();
93  }
94 
96 
97  #include "UEqn.H"
98  #include "TEqn.H"
99 
100  // --- Pressure corrector loop
101  while (pimple.correct())
102  {
103  #include "pEqn.H"
104  }
105 
106  if (pimple.turbCorr())
107  {
108  turbulence->correct();
109  }
110  }
111 
112  runTime.write();
113 
114  Info<< "ExecutionTime = "
115  << runTime.elapsedCpuTime()
116  << " s\n\n" << endl;
117  }
118 
119  Info<< "End\n" << endl;
120 
121  return 0;
122 }
123 
124 
125 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
interfaceProperties interface(alpha1, U, mixture())
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Read the control parameters used by setDeltaT.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
rhoEqn solve()
static const char nl
Definition: Ostream.H:262
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
rhoPhi
Definition: rhoEqn.H:10
messageStream Info
MULES: Multidimensional universal limiter for explicit solution.
Execute application functionObjects to post-process existing results.
Read the control parameters used by setDeltaT.