interPhaseChangeFoam.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  interPhaseChangeFoam
28 
29 Group
30  grpMultiphaseSolvers
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 
37  The momentum and other fluid properties are of the "mixture" and a
38  single momentum equation is solved.
39 
40  The set of phase-change models provided are designed to simulate cavitation
41  but other mechanisms of phase-change are supported within this solver
42  framework.
43 
44  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #include "fvCFD.H"
49 #include "CMULES.H"
50 #include "subCycle.H"
51 #include "interfaceProperties.H"
54 #include "pimpleControl.H"
55 #include "fvOptions.H"
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 int main(int argc, char *argv[])
60 {
61  argList::addNote
62  (
63  "Solver for two incompressible, isothermal immiscible fluids with"
64  " phase-change.\n"
65  "Uses VOF (volume of fluid) phase-fraction based interface capturing."
66  );
67 
68  #include "postProcess.H"
69 
70  #include "addCheckCaseOptions.H"
71  #include "setRootCaseLists.H"
72  #include "createTime.H"
73  #include "createMesh.H"
74  #include "createControl.H"
75  #include "createFields.H"
76  #include "createTimeControls.H"
77  #include "CourantNo.H"
78  #include "setInitialDeltaT.H"
79 
80  turbulence->validate();
81 
82  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83 
84  Info<< "\nStarting time loop\n" << endl;
85 
86  while (runTime.run())
87  {
88  #include "readTimeControls.H"
89  #include "CourantNo.H"
90  #include "setDeltaT.H"
91 
92  ++runTime;
93 
94  Info<< "Time = " << runTime.timeName() << nl << endl;
95 
96  // --- Pressure-velocity PIMPLE corrector loop
97  while (pimple.loop())
98  {
99  #include "alphaControls.H"
100 
102  (
103  IOobject
104  (
105  "rhoPhi",
106  runTime.timeName(),
107  mesh
108  ),
109  mesh,
111  );
112 
113  mixture->correct();
114 
115  #include "alphaEqnSubCycle.H"
116  interface.correct();
117 
118  #include "UEqn.H"
119 
120  // --- Pressure corrector loop
121  while (pimple.correct())
122  {
123  #include "pEqn.H"
124  }
125 
126  if (pimple.turbCorr())
127  {
128  turbulence->correct();
129  }
130  }
131 
132  runTime.write();
133 
134  runTime.printExecutionTime(Info);
135  }
136 
137  Info<< "End\n" << endl;
138 
139  return 0;
140 }
141 
142 
143 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
rhoPhi
rhoPhi
Definition: rhoEqn.H:10
CMULES.H
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
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
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
pimpleControl.H
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
addCheckCaseOptions.H
Required Classes.
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
interfaceProperties.H
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
postProcess.H
Execute application functionObjects to post-process existing results.
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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.
phaseChangeTwoPhaseMixture.H
readTimeControls.H
Read the control parameters used by setDeltaT.
createMesh.H
Required Variables.
createTime.H
fvCFD.H
mixture
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39