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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  interPhaseChangeDyMFoam
26 
27 Description
28  Solver for 2 incompressible, isothermal immiscible fluids with phase-change
29  (e.g. cavitation). Uses a VOF (volume of fluid) phase-fraction based
30  interface capturing approach, with optional mesh motion and mesh topology
31  changes including adaptive re-meshing.
32 
33  The momentum and other fluid properties are of the "mixture" and a
34  single momentum equation is solved.
35 
36  The set of phase-change models provided are designed to simulate cavitation
37  but other mechanisms of phase-change are supported within this solver
38  framework.
39 
40  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #include "fvCFD.H"
45 #include "dynamicFvMesh.H"
46 #include "CMULES.H"
47 #include "subCycle.H"
48 #include "interfaceProperties.H"
51 #include "pimpleControl.H"
52 #include "fvIOoptionList.H"
53 #include "CorrectPhi.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 int main(int argc, char *argv[])
59 {
60  #include "setRootCase.H"
61  #include "createTime.H"
62  #include "createDynamicFvMesh.H"
63 
64  pimpleControl pimple(mesh);
65 
66  #include "../interFoam/interDyMFoam/createControls.H"
67  #include "initContinuityErrs.H"
68  #include "createFields.H"
69  #include "createFvOptions.H"
70 
72  (
73  IOobject
74  (
75  "rAU",
76  runTime.timeName(),
77  mesh,
78  IOobject::READ_IF_PRESENT,
79  IOobject::AUTO_WRITE
80  ),
81  mesh,
82  dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
83  );
84 
85  #include "createUf.H"
86  #include "CourantNo.H"
87  #include "setInitialDeltaT.H"
88 
89  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
90 
91  Info<< "\nStarting time loop\n" << endl;
92 
93  while (runTime.run())
94  {
95  #include "../interFoam/interDyMFoam/readControls.H"
96 
97  // Store divU from the previous mesh so that it can be mapped
98  // and used in correctPhi to ensure the corrected phi has the
99  // same divergence
101 
102  #include "CourantNo.H"
103  #include "setDeltaT.H"
104 
105  runTime++;
106 
107  Info<< "Time = " << runTime.timeName() << nl << endl;
108 
109  // --- Pressure-velocity PIMPLE corrector loop
110  while (pimple.loop())
111  {
112  if (pimple.firstIter() || moveMeshOuterCorrectors)
113  {
114  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
115 
116  mesh.update();
117 
118  if (mesh.changing())
119  {
120  Info<< "Execution time for mesh.update() = "
121  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
122  << " s" << endl;
123 
124  gh = (g & mesh.C()) - ghRef;
125  ghf = (g & mesh.Cf()) - ghRef;
126  }
127 
128  if (mesh.changing() && correctPhi)
129  {
130  // Calculate absolute flux from the mapped surface velocity
131  phi = mesh.Sf() & Uf;
132 
133  #include "correctPhi.H"
134 
135  // Make the flux relative to the mesh motion
137  }
138 
139  if (mesh.changing() && checkMeshCourantNo)
140  {
141  #include "meshCourantNo.H"
142  }
143  }
144 
145  #include "alphaControls.H"
146 
148  (
149  IOobject
150  (
151  "rhoPhi",
152  runTime.timeName(),
153  mesh
154  ),
155  mesh,
157  );
158 
159  mixture->correct();
160 
161  #include "alphaEqnSubCycle.H"
162  interface.correct();
163 
164  #include "UEqn.H"
165 
166  // --- Pressure corrector loop
167  while (pimple.correct())
168  {
169  #include "pEqn.H"
170  }
171 
172  if (pimple.turbCorr())
173  {
174  turbulence->correct();
175  }
176  }
177 
178  runTime.write();
179 
180  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
181  << " ClockTime = " << runTime.elapsedClockTime() << " s"
182  << nl << endl;
183  }
184 
185  Info<< "End\n" << endl;
186 
187  return 0;
188 }
189 
190 
191 // ************************************************************************* //
Uf
Definition: pEqn.H:78
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Info<< "Reading field p_rgh\n"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating phaseChangeTwoPhaseMixture\n"<< endl;autoPtr< phaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:33
messageStream Info
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
dynamicFvMesh & mesh
rhoPhi
Definition: rhoEqn.H:10
const dictionary & pimple
checkMeshCourantNo
Definition: readControls.H:5
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:74
surfaceScalarField & phi
volScalarField rAU(1.0/UEqn.A())
const dimensionedVector & g
correctPhi
Definition: readControls.H:3
moveMeshOuterCorrectors
Definition: readControls.H:8
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Calculates and outputs the mean and maximum Courant Numbers.
int main(int argc, char *argv[])
Definition: postCalc.C:54
interfaceProperties interface(alpha1, U, mixture())
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
U
Definition: pEqn.H:82
const surfaceScalarField & ghf
Creates and initialises the velocity velocity field Uf.
const volScalarField & gh