compressibleInterFoam.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-2022 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  with optional mesh motion and mesh topology changes including adaptive
31  re-meshing.
32 
33  The momentum and other fluid properties are of the "mixture" and a single
34  momentum equation is solved.
35 
36  Either mixture or two-phase transport modelling may be selected. In the
37  mixture approach a single laminar, RAS or LES model is selected to model the
38  momentum stress. In the Euler-Euler two-phase approach separate laminar,
39  RAS or LES selected models are selected for each of the phases.
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #include "fvCFD.H"
44 #include "interfaceCompression.H"
45 #include "CMULES.H"
46 #include "EulerDdtScheme.H"
47 #include "localEulerDdtScheme.H"
48 #include "CrankNicolsonDdtScheme.H"
49 #include "subCycle.H"
51 #include "noPhaseChange.H"
52 #include "pimpleControl.H"
53 #include "pressureReference.H"
54 #include "fvModels.H"
55 #include "fvConstraints.H"
56 #include "CorrectPhi.H"
57 #include "fvcSmooth.H"
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 int main(int argc, char *argv[])
62 {
63  #include "postProcess.H"
64 
65  #include "setRootCaseLists.H"
66  #include "createTime.H"
67  #include "createMesh.H"
68  #include "initContinuityErrs.H"
69  #include "createDyMControls.H"
70  #include "createFields.H"
71  #include "createFieldRefs.H"
72  #include "CourantNo.H"
73  #include "setInitialDeltaT.H"
74  #include "createUfIfPresent.H"
75 
76  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77  Info<< "\nStarting time loop\n" << endl;
78 
79  while (pimple.run(runTime))
80  {
81  #include "readDyMControls.H"
82 
83  if (LTS)
84  {
85  #include "setRDeltaT.H"
86  }
87  else
88  {
89  #include "CourantNo.H"
90  #include "alphaCourantNo.H"
91  #include "setDeltaT.H"
92  }
93 
95 
96  // Store divU from the previous mesh so that it can be mapped
97  // and used in correctPhi to ensure the corrected phi has the
98  // same divergence
99  tmp<volScalarField> divU;
100 
101  if (correctPhi && mesh.topoChanged())
102  {
103  // Construct and register divU for mapping
104  divU = new volScalarField
105  (
106  "divU0",
108  );
109  }
110 
111  // Update the mesh for topology change, mesh to mesh mapping
112  mesh.update();
113 
114  runTime++;
115 
116  Info<< "Time = " << runTime.userTimeName() << nl << endl;
117 
118  // --- Pressure-velocity PIMPLE corrector loop
119  while (pimple.loop())
120  {
121  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
122  {
123  if (correctPhi && !divU.valid())
124  {
125  // Construct and register divU for correctPhi
126  divU = new volScalarField
127  (
128  "divU0",
130  );
131  }
132 
133  // Move the mesh
134  mesh.move();
135 
136  if (mesh.changing())
137  {
138  gh = (g & mesh.C()) - ghRef;
139  ghf = (g & mesh.Cf()) - ghRef;
140 
141  MRF.update();
142 
143  if (correctPhi)
144  {
145  #include "correctPhi.H"
146  }
147 
148  mixture.correct();
149 
150  if (checkMeshCourantNo)
151  {
152  #include "meshCourantNo.H"
153  }
154  }
155 
156  divU.clear();
157  }
158 
159  fvModels.correct();
160 
161  #include "alphaControls.H"
163 
164  turbulence.correctPhasePhi();
165 
166  #include "UEqn.H"
167  #include "TEqn.H"
168 
169  // --- Pressure corrector loop
170  while (pimple.correct())
171  {
172  #include "pEqn.H"
173  }
174 
175  if (pimple.turbCorr())
176  {
177  turbulence.correct();
178  }
179  }
180 
181  runTime.write();
182 
183  Info<< "ExecutionTime = "
184  << runTime.elapsedCpuTime()
185  << " s\n\n" << endl;
186  }
187 
188  Info<< "End\n" << endl;
189 
190  return 0;
191 }
192 
193 
194 // ************************************************************************* //
pimpleNoLoopControl & pimple
U
Definition: pEqn.H:72
IOMRFZoneList & MRF
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:353
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:255
correctPhi
checkMeshCourantNo
fvMesh & mesh
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
phi
Definition: correctPhi.H:3
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
const surfaceScalarField & ghf
static const char nl
Definition: Ostream.H:260
bool LTS
Definition: createRDeltaT.H:1
moveMeshOuterCorrectors
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
Foam::fvModels & fvModels
Calculates and outputs the mean and maximum Courant Numbers.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
const volScalarField & gh
messageStream Info
const dimensionedVector & g
Execute application functionObjects to post-process existing results.
Provides functions smooth spread and sweep which use the FvFaceCellWave algorithm to smooth and redis...
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Creates and initialises the velocity field Uf if required.