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  pimpleControl pimple(mesh);
86 
87  #include "createTimeControls.H"
88  #include "createDyMControls.H"
89  #include "initContinuityErrs.H"
90  #include "createFields.H"
91 
93  (
94  IOobject
95  (
96  "rAU",
97  runTime.timeName(),
98  mesh,
99  IOobject::READ_IF_PRESENT,
100  IOobject::AUTO_WRITE
101  ),
102  mesh,
103  dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
104  );
105 
106  #include "createUf.H"
107  #include "CourantNo.H"
108  #include "setInitialDeltaT.H"
109 
110  turbulence->validate();
111 
112  #include "setCellMask.H"
113  #include "setInterpolatedCells.H"
114 
115  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
116 
117  Info<< "\nStarting time loop\n" << endl;
118 
119  while (runTime.run())
120  {
121  #include "readControls.H"
122 
123  // Store divU from the previous mesh so that it can be mapped
124  // and used in correctPhi to ensure the corrected phi has the
125  // same divergence
127 
128  #include "CourantNo.H"
129  #include "setDeltaT.H"
130 
131  ++runTime;
132 
133  Info<< "Time = " << runTime.timeName() << nl << endl;
134 
135  // --- Pressure-velocity PIMPLE corrector loop
136  while (pimple.loop())
137  {
138  if (pimple.firstIter() || moveMeshOuterCorrectors)
139  {
140  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
141 
142  mesh.update();
143 
144  if (mesh.changing())
145  {
146  Info<< "Execution time for mesh.update() = "
147  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
148  << " s" << endl;
149 
150  gh = (g & mesh.C()) - ghRef;
151  ghf = (g & mesh.Cf()) - ghRef;
152 
153  // Update cellMask field for blocking out hole cells
154  #include "setCellMask.H"
155  #include "setInterpolatedCells.H"
156 
157  faceMask =
158  localMin<scalar>(mesh).interpolate(cellMask.oldTime());
159 
160 
161  // Zero Uf on old faceMask (H-I)
162  Uf *= faceMask;
163 
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  if (correctPhi)
178  {
179  #include "correctPhi.H"
180  }
181 
182  mixture->correct();
183 
184  // Zero phi on current H-I
185  faceMask = localMin<scalar>(mesh).interpolate(cellMask);
186 
187  phi *= faceMask;
188  U *= cellMask;
189 
190  // Make the flux relative to the mesh motion
192  }
193 
194  if (mesh.changing() && checkMeshCourantNo)
195  {
196  #include "meshCourantNo.H"
197  }
198  }
199 
200  #include "alphaControls.H"
201 
203  (
204  IOobject
205  (
206  "rhoPhi",
207  runTime.timeName(),
208  mesh
209  ),
210  mesh,
212  );
213 
214  mixture->correct();
215 
216  #include "alphaEqnSubCycle.H"
218  (
219  localMin<scalar>(mesh).interpolate(cellMask)
220  );
221  rhoPhi *= faceMask;
222 
223  interface.correct();
224 
225  #include "UEqn.H"
226 
227  // --- Pressure corrector loop
228  while (pimple.correct())
229  {
230  #include "pEqn.H"
231  }
232 
233  if (pimple.turbCorr())
234  {
235  turbulence->correct();
236  }
237  }
238 
239  runTime.write();
240 
241  runTime.printExecutionTime(Info);
242  }
243 
244  Info<< "End\n" << endl;
245 
246  return 0;
247 }
248 
249 
250 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
rhoPhi
rhoPhi
Definition: rhoEqn.H:10
CMULES.H
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
ghf
const surfaceScalarField & ghf
Definition: setRegionFluidFields.H:18
rAU
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
fvOptions.H
turbulence
Info<< "Reading field U\n"<< 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\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
interface
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
subCycle.H
turbulentTransportModel.H
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
faceMask
surfaceScalarField faceMask(localMin< scalar >(mesh).interpolate(cellMask))
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
createUf.H
Creates and initialises the velocity velocity field Uf.
setCellMask.H
Sets blocked cells mask field.
setInterpolatedCells.H
Sets blocked cells mask field.
correctPhi
correctPhi
Definition: readDyMControls.H:3
gh
const volScalarField & gh
Definition: setRegionFluidFields.H:17
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
rho
rho
Definition: readInitialConditions.H:88
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
interpolationCellPoint.H
pimpleControl.H
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:77
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
setRootCaseLists.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
interfaceProperties.H
divU
zeroField divU
Definition: alphaSuSp.H:3
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
postProcess.H
Execute application functionObjects to post-process existing results.
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
U
U
Definition: pEqn.H:72
localMin.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
createTimeControls.H
Read the control parameters used by setDeltaT.
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
phaseChangeTwoPhaseMixture.H
CorrectPhi.H
createTime.H
dynamicFvMesh.H
fvCFD.H
checkMeshCourantNo
checkMeshCourantNo
Definition: readDyMControls.H:9
transform.H
3D tensor transformation operations.
Foam::surfaceVectorField
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Definition: surfaceFieldsFwd.H:59
moveMeshOuterCorrectors
moveMeshOuterCorrectors
Definition: readDyMControls.H:15
cellCellStencilObject.H
mixture
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
oversetAdjustPhi.H
Adjust the balance of fluxes on the faces between interpolated and calculated to obey continuity.
createDynamicFvMesh.H
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190