interPhaseChangeDyMFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 "fvOptions.H"
53 #include "CorrectPhi.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 int main(int argc, char *argv[])
58 {
59  #include "postProcess.H"
60 
61  #include "setRootCaseLists.H"
62  #include "createTime.H"
63  #include "createDynamicFvMesh.H"
64  #include "createDyMControls.H"
65  #include "initContinuityErrs.H"
66  #include "createFields.H"
67 
69  (
70  IOobject
71  (
72  "rAU",
73  runTime.timeName(),
74  mesh,
75  IOobject::READ_IF_PRESENT,
76  IOobject::AUTO_WRITE
77  ),
78  mesh,
79  dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
80  );
81 
82  #include "createUf.H"
83  #include "CourantNo.H"
84  #include "setInitialDeltaT.H"
85 
86  turbulence->validate();
87 
88  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
89 
90  Info<< "\nStarting time loop\n" << endl;
91 
92  while (runTime.run())
93  {
94  #include "readDyMControls.H"
95 
96  // Store divU from the previous mesh so that it can be mapped
97  // and used in correctPhi to ensure the corrected phi has the
98  // same divergence
100 
101  #include "CourantNo.H"
102  #include "setDeltaT.H"
103 
104  runTime++;
105 
106  Info<< "Time = " << runTime.timeName() << nl << endl;
107 
108  // --- Pressure-velocity PIMPLE corrector loop
109  while (pimple.loop())
110  {
111  if (pimple.firstIter() || moveMeshOuterCorrectors)
112  {
113  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
114 
115  mesh.update();
116 
117  if (mesh.changing())
118  {
119  Info<< "Execution time for mesh.update() = "
120  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
121  << " s" << endl;
122 
123  gh = (g & mesh.C()) - ghRef;
124  ghf = (g & mesh.Cf()) - ghRef;
125  }
126 
127  if (mesh.changing() && correctPhi)
128  {
129  // Calculate absolute flux from the mapped surface velocity
130  phi = mesh.Sf() & Uf;
131 
132  #include "correctPhi.H"
133 
134  // Make the flux relative to the mesh motion
136  }
137 
138  if (mesh.changing() && checkMeshCourantNo)
139  {
140  #include "meshCourantNo.H"
141  }
142  }
143 
144  #include "alphaControls.H"
145 
147  (
148  IOobject
149  (
150  "rhoPhi",
151  runTime.timeName(),
152  mesh
153  ),
154  mesh,
156  );
157 
158  mixture->correct();
159 
160  #include "alphaEqnSubCycle.H"
161  interface.correct();
162 
163  #include "UEqn.H"
164 
165  // --- Pressure corrector loop
166  while (pimple.correct())
167  {
168  #include "pEqn.H"
169  }
170 
171  if (pimple.turbCorr())
172  {
173  turbulence->correct();
174  }
175  }
176 
177  runTime.write();
178 
179  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
180  << " ClockTime = " << runTime.elapsedClockTime() << " s"
181  << nl << endl;
182  }
183 
184  Info<< "End\n" << endl;
185 
186  return 0;
187 }
188 
189 
190 // ************************************************************************* //
pimpleNoLoopControl & pimple
surfaceScalarField & phi
interfaceProperties interface(alpha1, U, mixture())
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
correctPhi
checkMeshCourantNo
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
Creates and initialises the velocity velocity field Uf.
autoPtr< surfaceVectorField > Uf
const surfaceScalarField & ghf
static const char nl
Definition: Ostream.H:265
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
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating phaseChangeTwoPhaseMixture\"<< endl;autoPtr< phaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:33
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:188
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
rhoPhi
Definition: rhoEqn.H:10
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const volScalarField & gh
messageStream Info
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
const dimensionedVector & g
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.