overInterPhaseChangeDyMFoam.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 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
26 Application
27  overInterPhaseChangeDyMFoam
28 
29 Group
30  grpMultiphaseSolvers grpMovingMeshSolvers
31 
32 Description
33  Solver for two incompressible, isothermal, immiscible fluids with
34  phase-change (e.g. cavitation) using VOF (i.e. volume of fluid)
35  phase-fraction based interface capturing, with optional dynamic mesh
36  motion (including overset) and mesh topology changes including adaptive
37  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 #include "cellCellStencilObject.H"
62 #include "localMin.H"
63 #include "interpolationCellPoint.H"
64 #include "transform.H"
65 #include "oversetAdjustPhi.H"
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 int main(int argc, char *argv[])
70 {
71  argList::addNote
72  (
73  "Solver for two incompressible, isothermal, immiscible fluids with"
74  " phase-change\n"
75  "using VOF (volume of fluid) phase-fraction based interface capturing,"
76  " with optional dynamic mesh motion (including overset)\n"
77  "and mesh topology changes including adaptive re-meshing."
78  );
79 
80  #include "postProcess.H"
81 
82  #include "setRootCaseLists.H"
83  #include "createTime.H"
84  #include "createDynamicFvMesh.H"
85 
86  #include "createDyMControls.H"
87  #include "initContinuityErrs.H"
88  #include "createFields.H"
89 
91  (
92  IOobject
93  (
94  "rAU",
95  runTime.timeName(),
96  mesh,
97  IOobject::READ_IF_PRESENT,
98  IOobject::AUTO_WRITE
99  ),
100  mesh,
101  dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
102  );
103 
104  #include "createUf.H"
105  #include "CourantNo.H"
106  #include "setInitialDeltaT.H"
107 
108  turbulence->validate();
109 
110  #include "setCellMask.H"
111  #include "setInterpolatedCells.H"
112 
113  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
114 
115  Info<< "\nStarting time loop\n" << endl;
116 
117  while (runTime.run())
118  {
119  #include "readDyMControls.H"
120 
121  // Store divU from the previous mesh so that it can be mapped
122  // and used in correctPhi to ensure the corrected phi has the
123  // same divergence
125 
126  #include "CourantNo.H"
127  #include "setDeltaT.H"
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  gh = (g & mesh.C()) - ghRef;
149  ghf = (g & mesh.Cf()) - ghRef;
150 
151  // Update cellMask field for blocking out hole cells
152  #include "setCellMask.H"
153  #include "setInterpolatedCells.H"
154  #include "correctPhiFaceMask.H"
155 
156  mixture->correct();
157 
158  // Make the flux relative to the mesh motion
160  }
161 
162  if (mesh.changing() && checkMeshCourantNo)
163  {
164  #include "meshCourantNo.H"
165  }
166  }
167 
168  #include "alphaControls.H"
169 
171  (
172  IOobject
173  (
174  "rhoPhi",
175  runTime.timeName(),
176  mesh
177  ),
178  mesh,
180  );
181 
182  mixture->correct();
183 
184  #include "alphaEqnSubCycle.H"
185  rhoPhi *= faceMask;
186 
187  interface.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 // ************************************************************************* //
Sets blocked cells mask field.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
checkMeshCourantNo
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
Sets blocked cells mask field.
const surfaceScalarField & ghf
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
dynamicFvMesh & mesh
3D tensor transformation operations.
Creates and initialises the velocity velocity field Uf.
faceMask
Definition: setCellMask.H:43
const uniformDimensionedVectorField & g
pimpleControl & pimple
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
moveMeshOuterCorrectors
Correct phi on new faces C-I faces.
Calculates and outputs the mean and maximum Courant Numbers.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:183
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
rhoPhi
Definition: rhoEqn.H:10
const volScalarField & gh
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
zeroField divU
Definition: alphaSuSp.H:3
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:70
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Adjust the balance of fluxes on the faces between interpolated and calculated to obey continuity...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157