All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
compressibleInterFilmFoam.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-2020 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  compressibleInterFilmFoam
26 
27 Description
28  Solver for 2 compressible, non-isothermal immiscible fluids using a VOF
29  (volume of fluid) phase-fraction based interface capturing approach
30  and surface film modelling.
31 
32  The momentum and other fluid properties are of the "mixture" and a single
33  momentum equation is solved.
34 
35  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "fvCFD.H"
40 #include "CMULES.H"
41 #include "EulerDdtScheme.H"
42 #include "localEulerDdtScheme.H"
43 #include "CrankNicolsonDdtScheme.H"
44 #include "subCycle.H"
46 #include "pimpleControl.H"
47 #include "SLGThermo.H"
48 #include "surfaceFilmModel.H"
49 #include "fvOptions.H"
50 #include "fvcSmooth.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 "createTimeControls.H"
63  #include "createFields.H"
64  #include "createSurfaceFilmModel.H"
65 
66  volScalarField& p = mixture.p();
67  volScalarField& T = mixture.T();
68  const volScalarField& psi1 = mixture.thermo1().psi();
69  const volScalarField& psi2 = mixture.thermo2().psi();
70 
71  regionModels::surfaceFilmModel& surfaceFilm = tsurfaceFilm();
72 
73  if (!LTS)
74  {
75  #include "readTimeControls.H"
76  #include "CourantNo.H"
77  #include "setInitialDeltaT.H"
78  }
79 
80  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81 
82  Info<< "\nStarting time loop\n" << endl;
83 
84  while (pimple.run(runTime))
85  {
86  #include "readTimeControls.H"
87 
88  if (LTS)
89  {
90  #include "setRDeltaT.H"
91  }
92  else
93  {
94  #include "CourantNo.H"
95  #include "alphaCourantNo.H"
96  #include "setDeltaT.H"
97  }
98 
99  runTime++;
100 
101  Info<< "Time = " << runTime.timeName() << nl << endl;
102 
103  surfaceFilm.evolve();
104 
105  // --- Pressure-velocity PIMPLE corrector loop
106  while (pimple.loop())
107  {
108  #include "alphaControls.H"
110 
111  turbulence.correctPhasePhi();
112 
113  volScalarField::Internal Srho(surfaceFilm.Srho());
114  contErr -= posPart(Srho);
115 
116  #include "UEqn.H"
117  #include "TEqn.H"
118 
119  // --- Pressure corrector loop
120  while (pimple.correct())
121  {
122  #include "pEqn.H"
123  }
124 
125  if (pimple.turbCorr())
126  {
127  turbulence.correct();
128  }
129  }
130 
131  runTime.write();
132 
133  Info<< "ExecutionTime = "
134  << runTime.elapsedCpuTime()
135  << " s\n\n" << endl;
136  }
137 
138  Info<< "End\n" << endl;
139 
140  return 0;
141 }
142 
143 
144 // ************************************************************************* //
pimpleNoLoopControl & pimple
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Read the control parameters used by setDeltaT.
dimensionedScalar posPart(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
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
bool LTS
Definition: createRDeltaT.H:1
Info<< "\Constructing surface film model"<< endl;autoPtr< regionModels::surfaceFilmModel > tsurfaceFilm(regionModels::surfaceFilmModel::New(mesh, g))
messageStream Info
phaseChangeTwoPhaseMixture & mixture
Definition: createFields.H:38
Execute application functionObjects to post-process existing results.
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Read the control parameters used by setDeltaT.
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.