interPhaseChangeFoam.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  interPhaseChangeFoam
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.
31 
32  The momentum and other fluid properties are of the "mixture" and a
33  single momentum equation is solved.
34 
35  The set of phase-change models provided are designed to simulate cavitation
36  but other mechanisms of phase-change are supported within this solver
37  framework.
38 
39  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #include "fvCFD.H"
44 #include "CMULES.H"
45 #include "subCycle.H"
46 #include "interfaceProperties.H"
49 #include "pimpleControl.H"
50 #include "fvOptions.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 int main(int argc, char *argv[])
55 {
56  #include "postProcess.H"
57 
58  #include "setRootCaseLists.H"
59  #include "createTime.H"
60  #include "createMesh.H"
61  #include "createControl.H"
62  #include "createFields.H"
63  #include "createTimeControls.H"
64  #include "CourantNo.H"
65  #include "setInitialDeltaT.H"
66 
67  turbulence->validate();
68 
69  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71  Info<< "\nStarting time loop\n" << endl;
72 
73  while (runTime.run())
74  {
75  #include "readTimeControls.H"
76  #include "CourantNo.H"
77  #include "setDeltaT.H"
78 
79  runTime++;
80 
81  Info<< "Time = " << runTime.timeName() << nl << endl;
82 
83  // --- Pressure-velocity PIMPLE corrector loop
84  while (pimple.loop())
85  {
86  #include "alphaControls.H"
87 
89  (
90  IOobject
91  (
92  "rhoPhi",
93  runTime.timeName(),
94  mesh
95  ),
96  mesh,
98  );
99 
100  mixture->correct();
101 
102  #include "alphaEqnSubCycle.H"
103  interface.correct();
104 
105  #include "UEqn.H"
106 
107  // --- Pressure corrector loop
108  while (pimple.correct())
109  {
110  #include "pEqn.H"
111  }
112 
113  if (pimple.turbCorr())
114  {
115  turbulence->correct();
116  }
117  }
118 
119  runTime.write();
120 
121  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
122  << " ClockTime = " << runTime.elapsedClockTime() << " s"
123  << nl << endl;
124  }
125 
126  Info<< "End\n" << endl;
127 
128  return 0;
129 }
130 
131 
132 // ************************************************************************* //
pimpleNoLoopControl & pimple
interfaceProperties interface(alpha1, U, mixture())
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Read the control parameters used by setDeltaT.
dynamicFvMesh & mesh
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
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
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
messageStream Info
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Execute application functionObjects to post-process existing results.
Read the control parameters used by setDeltaT.
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.