interPhaseChangeDyMFoam.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 interPhaseChangeDyMFoam
28
29Group
30 grpMultiphaseSolvers grpMovingMeshSolvers
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 with optional mesh motion and mesh topology changes including
37 adaptive re-meshing.
38
39 The momentum and other fluid properties are of the "mixture" and a
40 single momentum equation is solved.
41
42 The set of phase-change models provided are designed to simulate cavitation
43 but other mechanisms of phase-change are supported within this solver
44 framework.
45
46 Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
47
48\*---------------------------------------------------------------------------*/
49
50#include "fvCFD.H"
51#include "dynamicFvMesh.H"
52#include "CMULES.H"
53#include "subCycle.H"
54#include "interfaceProperties.H"
57#include "pimpleControl.H"
58#include "fvOptions.H"
59#include "CorrectPhi.H"
60
61// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62
63int main(int argc, char *argv[])
64{
65 argList::addNote
66 (
67 "Solver for two incompressible, isothermal immiscible fluids with"
68 " phase-change.\n"
69 "Uses VOF (volume of fluid) phase-fraction based interface capturing,"
70 " with optional mesh motion and mesh topology changes including"
71 " adaptive re-meshing."
72 );
73
74 #include "postProcess.H"
75
76 #include "setRootCaseLists.H"
77 #include "createTime.H"
78 #include "createDynamicFvMesh.H"
79 #include "createDyMControls.H"
80 #include "initContinuityErrs.H"
81 #include "createFields.H"
82
84 (
85 IOobject
86 (
87 "rAU",
88 runTime.timeName(),
89 mesh,
90 IOobject::READ_IF_PRESENT,
91 IOobject::AUTO_WRITE
92 ),
93 mesh,
94 dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
95 );
96
97 #include "createUf.H"
98 #include "CourantNo.H"
99 #include "setInitialDeltaT.H"
100
101 turbulence->validate();
102
103 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104
105 Info<< "\nStarting time loop\n" << endl;
106
107 while (runTime.run())
108 {
109 #include "readDyMControls.H"
110
111 // Store divU from the previous mesh so that it can be mapped
112 // and used in correctPhi to ensure the corrected phi has the
113 // same divergence
114 volScalarField divU("divU0", fvc::div(fvc::absolute(phi, U)));
115
116 #include "CourantNo.H"
117 #include "setDeltaT.H"
118
119 ++runTime;
120
121 Info<< "Time = " << runTime.timeName() << nl << endl;
122
123 // --- Pressure-velocity PIMPLE corrector loop
124 while (pimple.loop())
125 {
126 if (pimple.firstIter() || moveMeshOuterCorrectors)
127 {
128 scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
129
130 mesh.update();
131
132 if (mesh.changing())
133 {
134 Info<< "Execution time for mesh.update() = "
135 << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
136 << " s" << endl;
137
138 gh = (g & mesh.C()) - ghRef;
139 ghf = (g & mesh.Cf()) - ghRef;
140 }
141
142 if (mesh.changing() && correctPhi)
143 {
144 // Calculate absolute flux from the mapped surface velocity
145 phi = mesh.Sf() & Uf;
146
147 #include "correctPhi.H"
148
149 // Make the flux relative to the mesh motion
150 fvc::makeRelative(phi, U);
151 }
152
153 if (mesh.changing() && checkMeshCourantNo)
154 {
155 #include "meshCourantNo.H"
156 }
157 }
158
159 #include "alphaControls.H"
160
162 (
163 IOobject
164 (
165 "rhoPhi",
166 runTime.timeName(),
167 mesh
168 ),
169 mesh,
170 dimensionedScalar(dimMass/dimTime, Zero)
171 );
172
173 mixture->correct();
174
175 #include "alphaEqnSubCycle.H"
176 interface.correct();
177
178 #include "UEqn.H"
179
180 // --- Pressure corrector loop
181 while (pimple.correct())
182 {
183 #include "pEqn.H"
184 }
185
186 if (pimple.turbCorr())
187 {
188 turbulence->correct();
189 }
190 }
191
192 runTime.write();
193
194 runTime.printExecutionTime(Info);
195 }
196
197 Info<< "End\n" << endl;
198
199 return 0;
200}
201
202
203// ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
rhoPhi
Definition: rhoEqn.H:10
const uniformDimensionedVectorField & g
surfaceScalarField & phi
const surfaceScalarField & ghf
const volScalarField & gh
pimpleControl & pimple
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
engineTime & runTime
autoPtr< surfaceVectorField > Uf
Creates and initialises the velocity velocity field Uf.
compressible::turbulenceModel & turbulence
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
zeroField divU
Definition: alphaSuSp.H:3
Calculates and outputs the mean and maximum Courant Numbers.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
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.
checkMeshCourantNo
correctPhi
moveMeshOuterCorrectors
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39