overCompressibleInterDyMFoam.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) 2021 OpenCFD Ltd.
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 overCompressibleInterDyMFoam
28
29Group
30 grpMultiphaseSolvers
31
32Description
33 Solver for two compressible, non-isothermal, immiscible fluids using VOF
34 (i.e. volume of fluid) phase-fraction based interface capturing approach.
35
36 This solver supports dynamic mesh motions including overset cases.
37
38 The momentum and other fluid properties are of the "mixture" and a single
39 momentum equation is solved.
40
41 Either mixture or two-phase transport modelling may be selected. In the
42 mixture approach, a single laminar, RAS or LES model is selected to model
43 the momentum stress. In the Euler-Euler two-phase approach separate
44 laminar, RAS or LES selected models are selected for each of the phases.
45
46\*---------------------------------------------------------------------------*/
47
48#include "fvCFD.H"
49#include "dynamicFvMesh.H"
50#include "CMULES.H"
51#include "EulerDdtScheme.H"
53#include "subCycle.H"
55#include "pimpleControl.H"
56#include "fvOptions.H"
57#include "fvcSmooth.H"
58
60#include "localMin.H"
62#include "transform.H"
63#include "fvMeshSubset.H"
64#include "oversetAdjustPhi.H"
65
66// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67
68int main(int argc, char *argv[])
69{
70 argList::addNote
71 (
72 "Solver for two compressible, non-isothermal, immiscible fluids"
73 " using VOF phase-fraction based interface capturing approach.\n"
74 "Supports dynamic mesh motions including overset cases."
75 );
76
77 #include "postProcess.H"
78
79 #include "addCheckCaseOptions.H"
80 #include "setRootCaseLists.H"
81 #include "createTime.H"
82 #include "createDynamicFvMesh.H"
83 pimpleControl pimple(mesh);
84
85 #include "createTimeControls.H"
86 #include "createDyMControls.H"
87 #include "createFields.H"
88
91 const volScalarField& psi1 = mixture.thermo1().psi();
92 const volScalarField& psi2 = mixture.thermo2().psi();
93
94 #include "correctPhi.H"
95 #include "createUf.H"
96
97 if (!LTS)
98 {
99 #include "CourantNo.H"
100 #include "setInitialDeltaT.H"
101 }
102
103 #include "setCellMask.H"
104 #include "setInterpolatedCells.H"
105
106 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
107
108 Info<< "\nStarting time loop\n" << endl;
109
110 while (runTime.run())
111 {
112 #include "readControls.H"
113
114 if (LTS)
115 {
116 #include "setRDeltaT.H"
117 }
118 else
119 {
120 #include "CourantNo.H"
121 #include "alphaCourantNo.H"
122 #include "setDeltaT.H"
123 }
124
125 ++runTime;
126
127 Info<< "Time = " << runTime.timeName() << nl << endl;
128
129 // --- Pressure-velocity PIMPLE corrector loop
130 while (pimple.loop())
131 {
132 if (pimple.firstIter() || moveMeshOuterCorrectors)
133 {
134 scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
135
136 mesh.update();
137
138 if (mesh.changing())
139 {
140 Info<< "Execution time for mesh.update() = "
141 << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
142 << " s" << endl;
143
144 // Do not apply previous time-step mesh compression flux
145 // if the mesh topology changed
146 if (mesh.topoChanging())
147 {
148 talphaPhi1Corr0.clear();
149 }
150
151 gh = (g & mesh.C()) - ghRef;
152 ghf = (g & mesh.Cf()) - ghRef;
153
154 // Update cellMask field for blocking out hole cells
155 #include "setCellMask.H"
156 #include "setInterpolatedCells.H"
157
158 faceMask =
159 localMin<scalar>(mesh).interpolate(cellMask.oldTime());
160
161 // Zero Uf on old faceMask (H-I)
162 Uf *= faceMask;
163
164 const surfaceVectorField Uint(fvc::interpolate(U));
165 // Update Uf and phi on new C-I faces
166 Uf += (1-faceMask)*Uint;
167
168 // Update Uf boundary
169 forAll(Uf.boundaryField(), patchI)
170 {
171 Uf.boundaryFieldRef()[patchI] =
172 Uint.boundaryField()[patchI];
173 }
174
175 phi = mesh.Sf() & Uf;
176
177 // Correct phi on individual regions
178 if (correctPhi)
179 {
180 #include "correctPhi.H"
181 }
182
183 mixture.correct();
184
185 // Zero phi on current H-I
186 faceMask = localMin<scalar>(mesh).interpolate(cellMask);
187
188 phi *= faceMask;
189 U *= cellMask;
190
191 // Make the flux relative to the mesh motion
192 fvc::makeRelative(phi, U);
193
194 }
195
196 if (mesh.changing() && checkMeshCourantNo)
197 {
198 #include "meshCourantNo.H"
199 }
200 }
201
202 #include "alphaControls.H"
203 #include "compressibleAlphaEqnSubCycle.H"
204
206 (
207 localMin<scalar>(mesh).interpolate(cellMask)
208 );
209 rhoPhi *= faceMask;
210
211 turbulence.correctPhasePhi();
212
213 #include "UEqn.H"
214 volScalarField divUp("divUp", fvc::div(fvc::absolute(phi, U), p));
215 #include "TEqn.H"
216
217 // --- Pressure corrector loop
218 while (pimple.correct())
219 {
220 #include "pEqn.H"
221 }
222
223 if (pimple.turbCorr())
224 {
225 turbulence.correct();
226 }
227 }
228
229 runTime.write();
230
231 runTime.printExecutionTime(Info);
232 }
233
234 Info<< "End\n" << endl;
235
236 return 0;
237}
238
239
240// ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Required Classes.
rhoPhi
Definition: rhoEqn.H:10
const uniformDimensionedVectorField & g
surfaceScalarField & phi
const surfaceScalarField & ghf
const volScalarField & gh
pimpleControl & pimple
const volScalarField & psi2
const volScalarField & psi1
U
Definition: pEqn.H:72
volScalarField & p
const volScalarField & T
surfaceScalarField faceMask(localMin< scalar >(mesh).interpolate(cellMask))
mesh interpolate(rAU)
tmp< surfaceScalarField > talphaPhi1Corr0
dynamicFvMesh & mesh
engineTime & runTime
bool LTS
Definition: createRDeltaT.H:1
Read the control parameters used by setDeltaT.
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...
Calculates and outputs the mean and maximum Courant Numbers.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
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
Adjust the balance of fluxes on the faces between interpolated and calculated to obey continuity.
Execute application functionObjects to post-process existing results.
checkMeshCourantNo
correctPhi
moveMeshOuterCorrectors
Sets blocked cells mask field.
Sets blocked cells mask field.
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
3D tensor transformation operations.