compressibleInterFilmFoam.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 | Copyright (C) 2011-2017 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  compressibleInterFoam
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 
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"
39 #include "CMULES.H"
40 #include "EulerDdtScheme.H"
41 #include "localEulerDdtScheme.H"
42 #include "CrankNicolsonDdtScheme.H"
43 #include "subCycle.H"
44 #include "rhoThermo.H"
45 #include "twoPhaseMixture.H"
46 #include "twoPhaseMixtureThermo.H"
48 #include "SLGThermo.H"
49 #include "surfaceFilmModel.H"
50 #include "pimpleControl.H"
51 #include "fvOptions.H"
52 #include "fvcSmooth.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 int main(int argc, char *argv[])
57 {
58  #include "postProcess.H"
59 
60  #include "setRootCase.H"
61  #include "createTime.H"
62  #include "createMesh.H"
63  #include "createControl.H"
64  #include "createTimeControls.H"
65  #include "createFields.H"
66  #include "createAlphaFluxes.H"
67  #include "createSurfaceFilmModel.H"
68  #include "createFvOptions.H"
69 
70  volScalarField& p = mixture.p();
71  volScalarField& T = mixture.T();
72  const volScalarField& psi1 = mixture.thermo1().psi();
73  const volScalarField& psi2 = mixture.thermo2().psi();
74 
75  filmModelType& surfaceFilm = tsurfaceFilm();
76 
77  turbulence->validate();
78 
79  if (!LTS)
80  {
81  #include "readTimeControls.H"
82  #include "CourantNo.H"
83  #include "setInitialDeltaT.H"
84  }
85 
86  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
87 
88  Info<< "\nStarting time loop\n" << endl;
89 
90  while (runTime.run())
91  {
92  #include "readTimeControls.H"
93 
94  if (LTS)
95  {
96  #include "setRDeltaT.H"
97  }
98  else
99  {
100  #include "CourantNo.H"
101  #include "alphaCourantNo.H"
102  #include "setDeltaT.H"
103  }
104 
105  runTime++;
106 
107  Info<< "Time = " << runTime.timeName() << nl << endl;
108 
109  surfaceFilm.evolve();
110 
111  // --- Pressure-velocity PIMPLE corrector loop
112  while (pimple.loop())
113  {
114  #include "alphaControls.H"
116 
117  volScalarField::Internal Srho(surfaceFilm.Srho());
118  contErr -= posPart(Srho);
119 
120  #include "UEqn.H"
121  #include "TEqn.H"
122 
123  // --- Pressure corrector loop
124  while (pimple.correct())
125  {
126  #include "pEqn.H"
127  }
128 
129  if (pimple.turbCorr())
130  {
131  turbulence->correct();
132  }
133  }
134 
135  runTime.write();
136 
137  Info<< "ExecutionTime = "
138  << runTime.elapsedCpuTime()
139  << " s\n\n" << endl;
140  }
141 
142  Info<< "End\n" << endl;
143 
144  return 0;
145 }
146 
147 
148 // ************************************************************************* //
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Read the control parameters used by setDeltaT.
dimensionedScalar posPart(const dimensionedScalar &ds)
Info<< "\Constructing surface film model"<< endl;typedef regionModels::surfaceFilmModels::surfaceFilmModel filmModelType;autoPtr< filmModelType > tsurfaceFilm(filmModelType::New(mesh, g))
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
static const char nl
Definition: Ostream.H:262
bool LTS
Definition: createRDeltaT.H:1
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
messageStream Info
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.