interPhaseChangeDyMFoam.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 -------------------------------------------------------------------------------
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  interPhaseChangeDyMFoam
28 
29 Group
30  grpMultiphaseSolvers grpMovingMeshSolvers
31 
32 Description
33  Solver for two incompressible, isothermal immiscible fluids with
34  phase-change (e.g. cavitation).
35  Uses VOF (volume of fluid) phase-fraction based interface capturing,
36  with optional mesh motion and mesh topology changes including
37  adaptive 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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 int main(int argc, char *argv[])
64 {
65  argList::addNote
66  (
67  "Solver for two incompressible, isothermal immiscible fluids with"
68  " phase-change.\n"
69  "Uses VOF (volume of fluid) phase-fraction based interface capturing,"
70  " with optional mesh motion and mesh topology changes including"
71  " adaptive re-meshing."
72  );
73 
74  #include "postProcess.H"
75 
76  #include "setRootCaseLists.H"
77  #include "createTime.H"
78  #include "createDynamicFvMesh.H"
79  #include "createDyMControls.H"
80  #include "initContinuityErrs.H"
81  #include "createFields.H"
82 
84  (
85  IOobject
86  (
87  "rAU",
88  runTime.timeName(),
89  mesh,
90  IOobject::READ_IF_PRESENT,
91  IOobject::AUTO_WRITE
92  ),
93  mesh,
94  dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
95  );
96 
97  #include "createUf.H"
98  #include "CourantNo.H"
99  #include "setInitialDeltaT.H"
100 
101  turbulence->validate();
102 
103  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104 
105  Info<< "\nStarting time loop\n" << endl;
106 
107  while (runTime.run())
108  {
109  #include "readDyMControls.H"
110 
111  // Store divU from the previous mesh so that it can be mapped
112  // and used in correctPhi to ensure the corrected phi has the
113  // same divergence
115 
116  #include "CourantNo.H"
117  #include "setDeltaT.H"
118 
119  ++runTime;
120 
121  Info<< "Time = " << runTime.timeName() << nl << endl;
122 
123  // --- Pressure-velocity PIMPLE corrector loop
124  while (pimple.loop())
125  {
126  if (pimple.firstIter() || moveMeshOuterCorrectors)
127  {
128  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
129 
130  mesh.update();
131 
132  if (mesh.changing())
133  {
134  Info<< "Execution time for mesh.update() = "
135  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
136  << " s" << endl;
137 
138  gh = (g & mesh.C()) - ghRef;
139  ghf = (g & mesh.Cf()) - ghRef;
140  }
141 
142  if (mesh.changing() && correctPhi)
143  {
144  // Calculate absolute flux from the mapped surface velocity
145  phi = mesh.Sf() & Uf;
146 
147  #include "correctPhi.H"
148 
149  // Make the flux relative to the mesh motion
151  }
152 
153  if (mesh.changing() && checkMeshCourantNo)
154  {
155  #include "meshCourantNo.H"
156  }
157  }
158 
159  #include "alphaControls.H"
160 
162  (
163  IOobject
164  (
165  "rhoPhi",
166  runTime.timeName(),
167  mesh
168  ),
169  mesh,
171  );
172 
173  mixture->correct();
174 
175  #include "alphaEqnSubCycle.H"
176  interface.correct();
177 
178  #include "UEqn.H"
179 
180  // --- Pressure corrector loop
181  while (pimple.correct())
182  {
183  #include "pEqn.H"
184  }
185 
186  if (pimple.turbCorr())
187  {
188  turbulence->correct();
189  }
190  }
191 
192  runTime.write();
193 
194  runTime.printExecutionTime(Info);
195  }
196 
197  Info<< "End\n" << endl;
198 
199  return 0;
200 }
201 
202 
203 // ************************************************************************* //
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
correctPhi
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
checkMeshCourantNo
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
const surfaceScalarField & ghf
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
dynamicFvMesh & mesh
Creates and initialises the velocity velocity field Uf.
autoPtr< surfaceVectorField > Uf
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
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.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157