compressibleInterDyMFoam.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 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Application
28 compressibleInterDyMFoam
29
30Description
31 Solver for two compressible, non-isothermal immiscible fluids using a VOF
32 (volume of fluid) phase-fraction based interface capturing approach,
33 with optional mesh motion and mesh topology changes including adaptive
34 re-meshing.
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"
44#include "dynamicFvMesh.H"
45#include "CMULES.H"
46#include "EulerDdtScheme.H"
47#include "localEulerDdtScheme.H"
49#include "subCycle.H"
51#include "pimpleControl.H"
52#include "fvOptions.H"
53#include "CorrectPhi.H"
54#include "fvcSmooth.H"
55
56// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57
58int main(int argc, char *argv[])
59{
60 argList::addNote
61 (
62 "Solver for two compressible, non-isothermal immiscible fluids"
63 " using VOF phase-fraction based interface capturing.\n"
64 "With optional mesh motion and mesh topology changes including"
65 " adaptive re-meshing."
66 );
67
68 #include "postProcess.H"
69
70 #include "setRootCaseLists.H"
71 #include "createTime.H"
72 #include "createDynamicFvMesh.H"
73 #include "initContinuityErrs.H"
74 #include "createDyMControls.H"
75 #include "createFields.H"
76 #include "createUf.H"
77 #include "CourantNo.H"
78 #include "setInitialDeltaT.H"
79
82 const volScalarField& psi1 = mixture.thermo1().psi();
83 const volScalarField& psi2 = mixture.thermo2().psi();
84
85 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 Info<< "\nStarting time loop\n" << endl;
87
88 while (runTime.run())
89 {
90 #include "readDyMControls.H"
91
92 // Store divU and divUp from the previous mesh so that it can be mapped
93 // and used in correctPhi to ensure the corrected phi has the
94 // same divergence
95 volScalarField divU("divU0", fvc::div(fvc::absolute(phi, U)));
96 volScalarField divUp("divUp", fvc::div(fvc::absolute(phi, U), p));
97
98 if (LTS)
99 {
100 #include "setRDeltaT.H"
101 }
102 else
103 {
104 #include "CourantNo.H"
105 #include "alphaCourantNo.H"
106 #include "setDeltaT.H"
107 }
108
109 ++runTime;
110
111 Info<< "Time = " << runTime.timeName() << nl << endl;
112
113 // --- Pressure-velocity PIMPLE corrector loop
114 while (pimple.loop())
115 {
116 if (pimple.firstIter() || moveMeshOuterCorrectors)
117 {
118 scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
119
120 mesh.update();
121
122 if (mesh.changing())
123 {
124 MRF.update();
125
126 Info<< "Execution time for mesh.update() = "
127 << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
128 << " s" << endl;
129
130 gh = (g & mesh.C()) - ghRef;
131 ghf = (g & mesh.Cf()) - ghRef;
132 }
133
134 if ((mesh.changing() && correctPhi))
135 {
136 // Calculate absolute flux from the mapped surface velocity
137 phi = mesh.Sf() & Uf;
138
139 #include "correctPhi.H"
140
141 // Make the fluxes relative to the mesh motion
142 fvc::makeRelative(phi, U);
143
144 mixture.correct();
145 }
146
147 if (mesh.changing() && checkMeshCourantNo)
148 {
149 #include "meshCourantNo.H"
150 }
151 }
152
153 #include "alphaControls.H"
154 #include "compressibleAlphaEqnSubCycle.H"
155
156 turbulence.correctPhasePhi();
157
158 #include "UEqn.H"
159 #include "TEqn.H"
160
161 // --- Pressure corrector loop
162 while (pimple.correct())
163 {
164 #include "pEqn.H"
165 }
166
167 if (pimple.turbCorr())
168 {
169 turbulence.correct();
170 }
171 }
172
174
175 // Correct p_rgh for consistency with p and the updated densities
176 p_rgh = p - rho*gh;
177 p_rgh.correctBoundaryConditions();
178
179 runTime.write();
180
181 runTime.printExecutionTime(Info);
182 }
183
184 Info<< "End\n" << endl;
185
186 return 0;
187}
188
189
190// ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
const uniformDimensionedVectorField & g
volScalarField & p_rgh
surfaceScalarField & phi
const surfaceScalarField & ghf
IOMRFZoneList & MRF
const volScalarField & gh
pimpleControl & pimple
const volScalarField & alpha1
const volScalarField & psi2
volScalarField & rho2
const volScalarField & alpha2
const volScalarField & psi1
volScalarField & rho1
U
Definition: pEqn.H:72
volScalarField & p
const volScalarField & T
dynamicFvMesh & mesh
engineTime & runTime
bool LTS
Definition: createRDeltaT.H:1
autoPtr< surfaceVectorField > Uf
Creates and initialises the velocity velocity field Uf.
compressible::turbulenceModel & turbulence
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
zeroField divU
Definition: alphaSuSp.H:3
Calculates and outputs the mean and maximum Courant Numbers.
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.
checkMeshCourantNo
correctPhi
moveMeshOuterCorrectors
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39