compressibleInterFilmFoam.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 compressibleInterFoam
28
29Description
30 Solver for two compressible, non-isothermal immiscible fluids using a VOF
31 (volume of fluid) phase-fraction based interface capturing approach.
32
33 The momentum and other fluid properties are of the "mixture" and a single
34 momentum equation is solved.
35
36 Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
37
38\*---------------------------------------------------------------------------*/
39
40#include "fvCFD.H"
41#include "CMULES.H"
42#include "EulerDdtScheme.H"
43#include "localEulerDdtScheme.H"
45#include "subCycle.H"
47#include "pimpleControl.H"
48#include "SLGThermo.H"
49#include "surfaceFilmModel.H"
50#include "pimpleControl.H"
51#include "fvOptions.H"
52#include "fvcSmooth.H"
53
54// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55
56int main(int argc, char *argv[])
57{
58 argList::addNote
59 (
60 "Solver for two compressible, non-isothermal immiscible fluids"
61 " using VOF phase-fraction based interface capturing."
62 );
63
64 #include "postProcess.H"
65
66 #include "addCheckCaseOptions.H"
67 #include "setRootCaseLists.H"
68 #include "createTime.H"
69 #include "createMesh.H"
70 #include "createControl.H"
71 #include "createTimeControls.H"
72 #include "createFields.H"
73 #include "createSurfaceFilmModel.H"
74
77 const volScalarField& psi1 = mixture.thermo1().psi();
78 const volScalarField& psi2 = mixture.thermo2().psi();
79
80 regionModels::surfaceFilmModel& surfaceFilm = tsurfaceFilm();
81
82 if (!LTS)
83 {
84 #include "readTimeControls.H"
85 #include "CourantNo.H"
86 #include "setInitialDeltaT.H"
87 }
88
89 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
90
91 Info<< "\nStarting time loop\n" << endl;
92
93 while (runTime.run())
94 {
95 #include "readTimeControls.H"
96
97 if (LTS)
98 {
99 #include "setRDeltaT.H"
100 }
101 else
102 {
103 #include "CourantNo.H"
104 #include "alphaCourantNo.H"
105 #include "setDeltaT.H"
106 }
107
108 ++runTime;
109
110 Info<< "Time = " << runTime.timeName() << nl << endl;
111
112 surfaceFilm.evolve();
113
114 // --- Pressure-velocity PIMPLE corrector loop
115 while (pimple.loop())
116 {
117 #include "alphaControls.H"
118 #include "compressibleAlphaEqnSubCycle.H"
119
120 turbulence.correctPhasePhi();
121
122 volScalarField::Internal Srho(surfaceFilm.Srho());
123 contErr -= posPart(Srho);
124
125 #include "UEqn.H"
126 #include "TEqn.H"
127
128 // --- Pressure corrector loop
129 while (pimple.correct())
130 {
131 #include "pEqn.H"
132 }
133
134 if (pimple.turbCorr())
135 {
136 turbulence.correct();
137 }
138 }
139
140 runTime.write();
141
142 runTime.printExecutionTime(Info);
143 }
144
145 Info<< "End\n" << endl;
146
147 return 0;
148}
149
150
151// ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Required Classes.
pimpleControl & pimple
const volScalarField & psi2
const volScalarField & psi1
volScalarField & p
const volScalarField & T
regionModels::surfaceFilmModel & surfaceFilm
Info<< "\nConstructing surface film model"<< endl;autoPtr< regionModels::surfaceFilmModel > tsurfaceFilm(regionModels::surfaceFilmModel::New(mesh, g))
volScalarField::Internal contErr((fvc::ddt(rho)+fvc::div(rhoPhi) -(fvOptions(alpha1, mixture.thermo1().rho())&rho1) -(fvOptions(alpha2, mixture.thermo2().rho())&rho2))())
engineTime & runTime
Required Variables.
bool LTS
Definition: createRDeltaT.H:1
Read the control parameters used by setDeltaT.
compressible::turbulenceModel & turbulence
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
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
dimensionedScalar posPart(const dimensionedScalar &ds)
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.
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39