compressibleMultiphaseInterFoam.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) 2013-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Application
27 compressibleMultiphaseInterFoam
28
29Group
30 grpMultiphaseSolvers
31
32Description
33 Solver for N compressible, non-isothermal immiscible fluids using a VOF
34 (volume of fluid) phase-fraction based interface capturing approach.
35
36 The momentum and other fluid properties are of the "mixture" and a single
37 momentum equation is solved.
38
39 Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
40
41\*---------------------------------------------------------------------------*/
42
43#include "fvCFD.H"
46#include "pimpleControl.H"
47
48// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49
50int main(int argc, char *argv[])
51{
52 argList::addNote
53 (
54 "Solver for N compressible, non-isothermal immiscible fluids"
55 " using VOF phase-fraction based interface capturing."
56 );
57
58 #include "postProcess.H"
59
60 #include "addCheckCaseOptions.H"
61 #include "setRootCaseLists.H"
62 #include "createTime.H"
63 #include "createMesh.H"
64 #include "createControl.H"
65 #include "createTimeControls.H"
66 #include "createFields.H"
67 #include "CourantNo.H"
68 #include "setInitialDeltaT.H"
69
72
73 turbulence->validate();
74
75 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76
77 Info<< "\nStarting time loop\n" << endl;
78
79 while (runTime.run())
80 {
81 #include "readTimeControls.H"
82 #include "CourantNo.H"
83 #include "alphaCourantNo.H"
84 #include "setDeltaT.H"
85
86 ++runTime;
87
88 Info<< "Time = " << runTime.timeName() << nl << endl;
89
90 // --- Pressure-velocity PIMPLE corrector loop
91 while (pimple.loop())
92 {
93 mixture.solve();
94
95 solve(fvm::ddt(rho) + fvc::div(mixture.rhoPhi()));
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 runTime.printExecutionTime(Info);
115 }
116
117 Info<< "End\n" << endl;
118
119 return 0;
120}
121
122
123// ************************************************************************* //
Required Classes.
pimpleControl & pimple
volScalarField & p
const volScalarField & T
engineTime & runTime
Required Variables.
Read the control parameters used by setDeltaT.
compressible::turbulenceModel & turbulence
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
Execute application functionObjects to post-process existing results.
Read the control parameters used by setDeltaT.
CEqn solve()
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39