compressibleInterFoam.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  Copyright (C) OpenCFD OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  compressibleInterFoam
29 
30 Group
31  grpMultiphaseSolvers
32 
33 Description
34  Solver for two compressible, non-isothermal immiscible fluids using a VOF
35  (volume of fluid) phase-fraction based interface capturing approach.
36 
37  The momentum and other fluid properties are of the "mixture" and a single
38  momentum equation is solved.
39 
40  Either mixture or two-phase transport modelling may be selected. In the
41  mixture approach a single laminar, RAS or LES model is selected to model the
42  momentum stress. In the Euler-Euler two-phase approach separate laminar,
43  RAS or LES selected models are selected for each of the phases.
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #include "fvCFD.H"
48 #include "CMULES.H"
49 #include "EulerDdtScheme.H"
50 #include "localEulerDdtScheme.H"
51 #include "CrankNicolsonDdtScheme.H"
52 #include "subCycle.H"
54 #include "pimpleControl.H"
55 #include "fvOptions.H"
56 #include "fvcSmooth.H"
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 int main(int argc, char *argv[])
61 {
62  argList::addNote
63  (
64  "Solver for two compressible, non-isothermal immiscible fluids"
65  " using VOF 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 "createTimeControls.H"
76  #include "createFields.H"
77 
78  volScalarField& p = mixture.p();
79  volScalarField& T = mixture.T();
80  const volScalarField& psi1 = mixture.thermo1().psi();
81  const volScalarField& psi2 = mixture.thermo2().psi();
82 
83  if (!LTS)
84  {
85  #include "readTimeControls.H"
86  #include "CourantNo.H"
87  #include "setInitialDeltaT.H"
88  }
89 
90  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91 
92  Info<< "\nStarting time loop\n" << endl;
93 
94  while (runTime.run())
95  {
96  #include "readTimeControls.H"
97 
98  if (LTS)
99  {
100  #include "setRDeltaT.H"
101  }
102  else
103  {
104  #include "CourantNo.H"
105  #include "alphaCourantNo.H"
106  #include "setDeltaT.H"
107  }
108 
109  ++runTime;
110 
111  Info<< "Time = " << runTime.timeName() << nl << endl;
112 
113  // --- Pressure-velocity PIMPLE corrector loop
114  while (pimple.loop())
115  {
116  #include "alphaControls.H"
117  #include "compressibleAlphaEqnSubCycle.H"
118 
119  turbulence.correctPhasePhi();
120 
121  #include "UEqn.H"
122  volScalarField divUp("divUp", fvc::div(fvc::absolute(phi, U), p));
123  #include "TEqn.H"
124 
125  // --- Pressure corrector loop
126  while (pimple.correct())
127  {
128  #include "pEqn.H"
129  }
130 
131  if (pimple.turbCorr())
132  {
133  turbulence.correct();
134  }
135  }
136 
137  runTime.write();
138 
139  runTime.printExecutionTime(Info);
140  }
141 
142  Info<< "End\n" << endl;
143 
144  return 0;
145 }
146 
147 
148 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
p
volScalarField & p
Definition: createFieldRefs.H:8
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
subCycle.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
psi2
const volScalarField & psi2
Definition: setRegionFluidFields.H:31
pimpleControl.H
compressibleInterPhaseTransportModel.H
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
addCheckCaseOptions.H
Required Classes.
localEulerDdtScheme.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
EulerDdtScheme.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
CrankNicolsonDdtScheme.H
LTS
bool LTS
Definition: createRDeltaT.H:1
postProcess.H
Execute application functionObjects to post-process existing results.
psi1
const volScalarField & psi1
Definition: setRegionFluidFields.H:28
T
const volScalarField & T
Definition: createFieldRefs.H:2
U
U
Definition: pEqn.H:72
Foam::nl
constexpr char nl
Definition: Ostream.H:404
createTimeControls.H
Read the control parameters used by setDeltaT.
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
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
fvcSmooth.H
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...