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 Group
28  grpMultiphaseSolvers grpMovingMeshSolvers
29 
30 Description
31  Solver for 2 incompressible, isothermal immiscible fluids with phase-change
32  (e.g. cavitation). Uses a VOF (volume of fluid) phase-fraction based
33  interface capturing approach, with optional mesh motion and mesh topology
34  changes including adaptive re-meshing.
35 
36  The momentum and other fluid properties are of the "mixture" and a
37  single momentum equation is solved.
38 
39  The set of phase-change models provided are designed to simulate cavitation
40  but other mechanisms of phase-change are supported within this solver
41  framework.
42 
43  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #include "fvCFD.H"
48 #include "dynamicFvMesh.H"
49 #include "CMULES.H"
50 #include "subCycle.H"
51 #include "interfaceProperties.H"
54 #include "pimpleControl.H"
55 #include "fvOptions.H"
56 #include "CorrectPhi.H"
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 int main(int argc, char *argv[])
62 {
63  #include "setRootCase.H"
64  #include "createTime.H"
65  #include "createDynamicFvMesh.H"
66 
67  pimpleControl pimple(mesh);
68 
69  #include "../interFoam/interDyMFoam/createControls.H"
70  #include "initContinuityErrs.H"
71  #include "createFields.H"
72  #include "createFvOptions.H"
73 
75  (
76  IOobject
77  (
78  "rAU",
79  runTime.timeName(),
80  mesh,
81  IOobject::READ_IF_PRESENT,
82  IOobject::AUTO_WRITE
83  ),
84  mesh,
85  dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
86  );
87 
88  #include "createUf.H"
89  #include "CourantNo.H"
90  #include "setInitialDeltaT.H"
91 
92  turbulence->validate();
93 
94  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
95 
96  Info<< "\nStarting time loop\n" << endl;
97 
98  while (runTime.run())
99  {
100  #include "../interFoam/interDyMFoam/readControls.H"
101 
102  // Store divU from the previous mesh so that it can be mapped
103  // and used in correctPhi to ensure the corrected phi has the
104  // same divergence
106 
107  #include "CourantNo.H"
108  #include "setDeltaT.H"
109 
110  runTime++;
111 
112  Info<< "Time = " << runTime.timeName() << nl << endl;
113 
114  // --- Pressure-velocity PIMPLE corrector loop
115  while (pimple.loop())
116  {
117  if (pimple.firstIter() || moveMeshOuterCorrectors)
118  {
119  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
120 
121  mesh.update();
122 
123  if (mesh.changing())
124  {
125  Info<< "Execution time for mesh.update() = "
126  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
127  << " s" << endl;
128 
129  gh = (g & mesh.C()) - ghRef;
130  ghf = (g & mesh.Cf()) - ghRef;
131  }
132 
133  if (mesh.changing() && correctPhi)
134  {
135  // Calculate absolute flux from the mapped surface velocity
136  phi = mesh.Sf() & Uf;
137 
138  #include "correctPhi.H"
139 
140  // Make the flux relative to the mesh motion
142  }
143 
144  if (mesh.changing() && checkMeshCourantNo)
145  {
146  #include "meshCourantNo.H"
147  }
148  }
149 
150  #include "alphaControls.H"
151 
153  (
154  IOobject
155  (
156  "rhoPhi",
157  runTime.timeName(),
158  mesh
159  ),
160  mesh,
162  );
163 
164  mixture->correct();
165 
166  #include "alphaEqnSubCycle.H"
167  interface.correct();
168 
169  #include "UEqn.H"
170 
171  // --- Pressure corrector loop
172  while (pimple.correct())
173  {
174  #include "pEqn.H"
175  }
176 
177  if (pimple.turbCorr())
178  {
179  turbulence->correct();
180  }
181  }
182 
183  runTime.write();
184 
185  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
186  << " ClockTime = " << runTime.elapsedClockTime() << " s"
187  << nl << endl;
188  }
189 
190  Info<< "End\n" << endl;
191 
192  return 0;
193 }
194 
195 
196 // ************************************************************************* //
rhoPhi
rhoPhi
Definition: rhoEqn.H:10
CMULES.H
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
fvOptions.H
ghf
const surfaceScalarField & ghf
Definition: setRegionFluidFields.H:35
subCycle.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
turbulentTransportModel.H
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
moveMeshOuterCorrectors
moveMeshOuterCorrectors
Definition: readControls.H:8
createUf.H
Creates and initialises the velocity velocity field Uf.
mixture
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
pimpleControl.H
U
U
Definition: pEqn.H:46
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:74
createFvOptions.H
interface
interfaceProperties interface(alpha1, U, mixture())
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
checkMeshCourantNo
bool checkMeshCourantNo
Definition: readControls.H:6
interfaceProperties.H
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
gh
const volScalarField & gh
Definition: setRegionFluidFields.H:34
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
rho
rho
Definition: pEqn.H:3
rAU
volScalarField rAU("rAU", 1.0/UEqn().A())
setRootCase.H
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
phaseChangeTwoPhaseMixture.H
CorrectPhi.H
createTime.H
dynamicFvMesh.H
fvCFD.H
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
fixedFluxPressureFvPatchScalarField.H
Uf
Uf
Definition: pEqn.H:78
createDynamicFvMesh.H
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
correctPhi
bool correctPhi
Definition: readControls.H:3