compressibleMultiphaseInterFoam.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) 2013-2021 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  compressibleMultiphaseInterFoam
26 
27 Description
28  Solver for n compressible, non-isothermal immiscible fluids using a VOF
29  (volume of fluid) phase-fraction based interface capturing approach.
30 
31  The momentum and other fluid properties are of the "mixture" and a single
32  momentum equation is solved.
33 
34  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "fvCFD.H"
41 #include "pimpleControl.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 int main(int argc, char *argv[])
46 {
47  #include "postProcess.H"
48 
49  #include "setRootCaseLists.H"
50  #include "createTime.H"
51  #include "createMesh.H"
52  #include "createControl.H"
53  #include "createTimeControls.H"
54  #include "createFields.H"
55  #include "CourantNo.H"
56  #include "setInitialDeltaT.H"
57 
58  volScalarField& p = mixture.p();
59  volScalarField& T = mixture.T();
60 
61  turbulence->validate();
62 
63  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65  Info<< "\nStarting time loop\n" << endl;
66 
67  while (pimple.run(runTime))
68  {
69  #include "readTimeControls.H"
70  #include "CourantNo.H"
71  #include "alphaCourantNo.H"
72  #include "setDeltaT.H"
73 
74  runTime++;
75 
76  Info<< "Time = " << runTime.timeName() << nl << endl;
77 
78  // --- Pressure-velocity PIMPLE corrector loop
79  while (pimple.loop())
80  {
81  mixture.solve();
82 
83  solve(fvm::ddt(rho) + fvc::div(mixture.rhoPhi()));
84 
85  #include "UEqn.H"
86  #include "TEqn.H"
87 
88  // --- Pressure corrector loop
89  while (pimple.correct())
90  {
91  #include "pEqn.H"
92  }
93 
94  if (pimple.turbCorr())
95  {
96  turbulence->correct();
97  }
98  }
99 
100  runTime.write();
101 
102  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
103  << " ClockTime = " << runTime.elapsedClockTime() << " s"
104  << nl << endl;
105  }
106 
107  Info<< "End\n" << endl;
108 
109  return 0;
110 }
111 
112 
113 // ************************************************************************* //
pimpleNoLoopControl & pimple
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:251
Read the control parameters used by setDeltaT.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
static const char nl
Definition: Ostream.H:260
rhoEqn solve()
messageStream Info
Execute application functionObjects to post-process existing results.
Read the control parameters used by setDeltaT.