compressibleInterIsoFoam.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) 2020 OpenCFD Ltd.
9  Copyright (C) 2020 Johan Roenby
10  Copyright (C) 2020 DLR
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 Application
29  compressibleInterFlow
30 
31 Description
32  Solver derived from interFoam for two compressible, immiscible
33  fluids using the isoAdvector phase-fraction based interface capturing
34  approach, with optional mesh motion and mesh topology changes including
35  adaptive re-meshing.
36 
37  Reference:
38  \verbatim
39  Roenby, J., Bredmose, H. and Jasak, H. (2016).
40  A computational method for sharp interface advection
41  Royal Society Open Science, 3
42  doi 10.1098/rsos.160405
43  \endverbatim
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #include "fvCFD.H"
48 #include "dynamicFvMesh.H"
49 #include "CMULES.H"
50 #include "EulerDdtScheme.H"
51 #include "localEulerDdtScheme.H"
52 #include "CrankNicolsonDdtScheme.H"
53 #include "subCycle.H"
55 #include "pimpleControl.H"
56 #include "fvOptions.H"
57 #include "CorrectPhi.H"
58 #include "fvcSmooth.H"
59 #include "dynamicRefineFvMesh.H"
60 #include "isoAdvection.H"
61 #include "twoPhaseMixtureThermo.H"
62 
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 int main(int argc, char *argv[])
67 {
68  argList::addNote
69  (
70  "Solver for two compressible, non-isothermal immiscible fluids"
71  " using VOF phase-fraction based interface capturing.\n"
72  "With optional mesh motion and mesh topology changes including"
73  " adaptive re-meshing."
74  );
75 
76  #include "postProcess.H"
77 
78  #include "setRootCaseLists.H"
79  #include "createTime.H"
80  #include "createDynamicFvMesh.H"
81  #include "initContinuityErrs.H"
82  #include "createDyMControls.H"
83  #include "createFields.H"
84  #include "createUf.H"
85  #include "CourantNo.H"
86  #include "setInitialDeltaT.H"
87 
88  volScalarField& p = mixture.p();
89  volScalarField& T = mixture.T();
90  const volScalarField& psi1 = mixture.thermo1().psi();
91  const volScalarField& psi2 = mixture.thermo2().psi();
92 
93  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
94  Info<< "\nStarting time loop\n" << endl;
95 
96  while (runTime.run())
97  {
98  #include "readDyMControls.H"
99 
100  // Store divU and divUp from the previous mesh so that it can be mapped
101  // and used in correctPhi to ensure the corrected phi has the
102  // same divergence
104 
105  #include "CourantNo.H"
106  #include "alphaCourantNo.H"
107  #include "setDeltaT.H"
108 
109 
110  ++runTime;
111 
112  Info<< "Time = " << runTime.timeName() << nl << endl;
113 
114  // --- Pressure-velocity PIMPLE corrector loop
115  while (pimple.loop())
116  {
117  if (pimple.firstIter() || moveMeshOuterCorrectors)
118  {
119  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
120 
121  if (isA<dynamicRefineFvMesh>(mesh))
122  {
123  advector.surf().reconstruct();
124  }
125 
126  mesh.update();
127 
128  if (mesh.changing())
129  {
130  gh = (g & mesh.C()) - ghRef;
131  ghf = (g & mesh.Cf()) - ghRef;
132 
133  if (isA<dynamicRefineFvMesh>(mesh))
134  {
135  advector.surf().mapAlphaField();
136  alpha2 = 1.0 - alpha1;
137  alpha2.correctBoundaryConditions();
138  rho == alpha1*rho1 + alpha2*rho2;
139  rho.correctBoundaryConditions();
140  rho.oldTime() = rho;
141  alpha2.oldTime() = alpha2;
142  }
143 
144  MRF.update();
145 
146  Info<< "Execution time for mesh.update() = "
147  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
148  << " s" << endl;
149 
150  }
151 
152  if ((mesh.changing() && correctPhi))
153  {
154  // Calculate absolute flux from the mapped surface velocity
155  phi = mesh.Sf() & Uf;
156 
157  #include "correctPhi.H"
158 
159  // Make the fluxes relative to the mesh motion
161 
162  mixture.correct();
163  }
164 
165  if (mesh.changing() && checkMeshCourantNo)
166  {
167  #include "meshCourantNo.H"
168  }
169  }
170 
171  #include "alphaControls.H"
172  #include "compressibleAlphaEqnSubCycle.H"
173 
174  turbulence.correctPhasePhi();
175 
176  #include "UEqn.H"
177  volScalarField divUp("divUp", fvc::div(fvc::absolute(phi, U), p));
178  #include "TEqn.H"
179 
180  // --- Pressure corrector loop
181  while (pimple.correct())
182  {
183  #include "pEqn.H"
184  }
185 
186  if (pimple.turbCorr())
187  {
188  turbulence.correct();
189  }
190  }
191 
192  runTime.write();
193 
194  runTime.printExecutionTime(Info);
195  }
196 
197  Info<< "End\n" << endl;
198 
199  return 0;
200 }
201 
202 
203 // ************************************************************************* //
isoAdvection advector(alpha1, phi, U)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
volScalarField & rho1
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
correctPhi
const volScalarField & alpha2
checkMeshCourantNo
const surfaceScalarField & ghf
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
dynamicFvMesh & mesh
const volScalarField & psi2
IOMRFZoneList & MRF
Creates and initialises the velocity velocity field Uf.
autoPtr< surfaceVectorField > Uf
const uniformDimensionedVectorField & g
pimpleControl & pimple
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
const volScalarField & psi1
volScalarField & rho2
moveMeshOuterCorrectors
const volScalarField & T
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:183
U
Definition: pEqn.H:72
const volScalarField & gh
messageStream Info
Information stream (stdout output on master, null elsewhere)
zeroField divU
Definition: alphaSuSp.H:3
volScalarField & p
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:70
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
const volScalarField & alpha1