interPhaseChangeFoam.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-------------------------------------------------------------------------------
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 interPhaseChangeFoam
28
29Group
30 grpMultiphaseSolvers
31
32Description
33 Solver for two incompressible, isothermal immiscible fluids with
34 phase-change (e.g. cavitation).
35 Uses VOF (volume of fluid) phase-fraction based interface capturing.
36
37 The momentum and other fluid properties are of the "mixture" and a
38 single momentum equation is solved.
39
40 The set of phase-change models provided are designed to simulate cavitation
41 but other mechanisms of phase-change are supported within this solver
42 framework.
43
44 Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
45
46\*---------------------------------------------------------------------------*/
47
48#include "fvCFD.H"
49#include "CMULES.H"
50#include "subCycle.H"
51#include "interfaceProperties.H"
54#include "pimpleControl.H"
55#include "fvOptions.H"
56
57// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58
59int main(int argc, char *argv[])
60{
61 argList::addNote
62 (
63 "Solver for two incompressible, isothermal immiscible fluids with"
64 " phase-change.\n"
65 "Uses VOF (volume of fluid) phase-fraction based interface capturing."
66 );
67
68 #include "postProcess.H"
69
70 #include "addCheckCaseOptions.H"
71 #include "setRootCaseLists.H"
72 #include "createTime.H"
73 #include "createMesh.H"
74 #include "createControl.H"
75 #include "createFields.H"
76 #include "createTimeControls.H"
77 #include "CourantNo.H"
78 #include "setInitialDeltaT.H"
79
80 turbulence->validate();
81
82 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83
84 Info<< "\nStarting time loop\n" << endl;
85
86 while (runTime.run())
87 {
88 #include "readTimeControls.H"
89 #include "CourantNo.H"
90 #include "setDeltaT.H"
91
92 ++runTime;
93
94 Info<< "Time = " << runTime.timeName() << nl << endl;
95
96 // --- Pressure-velocity PIMPLE corrector loop
97 while (pimple.loop())
98 {
99 #include "alphaControls.H"
100
102 (
103 IOobject
104 (
105 "rhoPhi",
106 runTime.timeName(),
107 mesh
108 ),
109 mesh,
110 dimensionedScalar(dimMass/dimTime, Zero)
111 );
112
113 mixture->correct();
114
115 #include "alphaEqnSubCycle.H"
116 interface.correct();
117
118 #include "UEqn.H"
119
120 // --- Pressure corrector loop
121 while (pimple.correct())
122 {
123 #include "pEqn.H"
124 }
125
126 if (pimple.turbCorr())
127 {
128 turbulence->correct();
129 }
130 }
131
132 runTime.write();
133
134 runTime.printExecutionTime(Info);
135 }
136
137 Info<< "End\n" << endl;
138
139 return 0;
140}
141
142
143// ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Required Classes.
rhoPhi
Definition: rhoEqn.H:10
pimpleControl & pimple
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
Read the control parameters used by setDeltaT.
compressible::turbulenceModel & turbulence
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
Information stream (stdout output on master, null elsewhere)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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.
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39