Loading...
Searching...
No Matches
overInterDyMFoam.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-2014 OpenFOAM Foundation
9 Copyright (C) 2016-2021 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 overInterDyMFoam
29
30Group
31 grpMultiphaseSolvers grpMovingMeshSolvers
32
33Description
34 Solver for two incompressible, isothermal immiscible fluids using a VOF
35 (volume of fluid) phase-fraction based interface capturing approach,
36 with optional mesh motion and mesh topology changes including adaptive
37 re-meshing.
38
39\*---------------------------------------------------------------------------*/
40
41#include "fvCFD.H"
42#include "dynamicFvMesh.H"
43#include "CMULES.H"
44#include "EulerDdtScheme.H"
45#include "localEulerDdtScheme.H"
47#include "subCycle.H"
50#include "pimpleControl.H"
51#include "fvOptions.H"
52#include "fvcSmooth.H"
54#include "localMin.H"
55#include "oversetAdjustPhi.H"
56#include "oversetPatchPhiErr.H"
57
58// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59
60int main(int argc, char *argv[])
61{
62 argList::addNote
63 (
64 "Solver for two incompressible, isothermal immiscible fluids using"
65 " VOF phase-fraction based interface capturing\n"
66 "With optional mesh motion and mesh topology changes including"
67 " adaptive re-meshing."
68 );
69
70 #include "postProcess.H"
71
72 #include "setRootCaseLists.H"
73 #include "createTime.H"
74 #include "createDynamicFvMesh.H"
75 #include "initContinuityErrs.H"
76
77 #include "createDyMControls.H"
78 #include "createFields.H"
79 #include "createAlphaFluxes.H"
80 #include "createFvOptions.H"
81
83 (
84 IOobject
85 (
86 "rAU",
87 runTime.timeName(),
88 mesh,
89 IOobject::READ_IF_PRESENT,
90 IOobject::AUTO_WRITE
91 ),
92 mesh,
93 dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
94 );
95
96 #include "createUf.H"
97 #include "createControls.H"
98
99 #include "setCellMask.H"
100 #include "setInterpolatedCells.H"
101
102 turbulence->validate();
103
104 if (!LTS)
105 {
106 #include "CourantNo.H"
107 #include "setInitialDeltaT.H"
108 }
109
110 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
111 Info<< "\nStarting time loop\n" << endl;
112
113 while (runTime.run())
114 {
115 #include "readDyMControls.H"
116 #include "readOversetDyMControls.H"
117
118 if (LTS)
119 {
120 #include "setRDeltaT.H"
121 }
122 else
123 {
124 #include "CourantNo.H"
125 #include "alphaCourantNo.H"
126 #include "setDeltaT.H"
127 }
128
129 ++runTime;
130
131 Info<< "Time = " << runTime.timeName() << nl << endl;
132
133 // --- Pressure-velocity PIMPLE corrector loop
134 while (pimple.loop())
135 {
136 if (pimple.firstIter() || moveMeshOuterCorrectors)
137 {
138 scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
139
140 mesh.update();
141
142 if (mesh.changing())
143 {
144 Info<< "Execution time for mesh.update() = "
145 << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
146 << " s" << endl;
147
148 // Do not apply previous time-step mesh compression flux
149 // if the mesh topology changed
150 if (mesh.topoChanging())
151 {
152 talphaPhi1Corr0.clear();
153 }
154
155 // Update cellMask field for blocking out hole cells
156 #include "setCellMask.H"
157 #include "setInterpolatedCells.H"
158 #include "correctPhiFaceMask.H"
159
160 gh = (g & mesh.C()) - ghRef;
161 ghf = (g & mesh.Cf()) - ghRef;
162
163
164 mixture.correct();
165
166 // Make the flux relative to the mesh motion
167 fvc::makeRelative(phi, U);
168
169 }
170
171 if (mesh.changing() && checkMeshCourantNo)
172 {
173 #include "meshCourantNo.H"
174 }
175 }
176
177 if (adjustFringe)
178 {
180 }
181
182 #include "alphaControls.H"
183 #include "alphaEqnSubCycle.H"
184
185 rhoPhi *= faceMask;
186
187 mixture.correct();
188
189 #include "UEqn.H"
190
191 // --- Pressure corrector loop
192 while (pimple.correct())
193 {
194 #include "pEqn.H"
195 }
196
197 if (pimple.turbCorr())
198 {
199 turbulence->correct();
200 }
201 }
202
203 runTime.write();
204
205 runTime.printExecutionTime(Info);
206 }
207
208 Info<< "End\n" << endl;
209
210 return 0;
211}
212
213
214// ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
rhoPhi
Definition rhoEqn.H:10
Calculates and outputs the mean and maximum Courant Numbers.
const uniformDimensionedVectorField & g
const surfaceScalarField & ghf
const volScalarField & gh
pimpleControl & pimple
U
Definition pEqn.H:72
Correct phi on new faces C-I faces.
tmp< surfaceScalarField > talphaPhi1Corr0
dynamicFvMesh & mesh
engineTime & runTime
bool LTS
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...
tmp< volScalarField > rAU
Calculates and outputs the mean and maximum Courant Numbers.
Calculates and outputs the mean and maximum Courant Numbers.
bool oversetAdjustPhi(surfaceScalarField &phi, const volVectorField &U, const label zoneId=-1)
Adjust the balance of fluxes to obey continuity.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Adjust the balance of fluxes on the faces between interpolated and calculated to obey continuity.
oversetPatchPhiErr
Execute application functionObjects to post-process existing results.
checkMeshCourantNo
moveMeshOuterCorrectors
label zoneIdMass
Sets blocked cells mask field.
faceMask
Definition setCellMask.H:43
Sets blocked cells mask field.
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture