compressibleInterDyMFoam.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-2018 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  compressibleInterDyMFoam
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  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #include "fvCFD.H"
41 #include "dynamicFvMesh.H"
42 #include "CMULES.H"
43 #include "EulerDdtScheme.H"
44 #include "localEulerDdtScheme.H"
45 #include "CrankNicolsonDdtScheme.H"
46 #include "subCycle.H"
48 #include "pimpleControl.H"
49 #include "fvOptions.H"
50 #include "CorrectPhi.H"
51 #include "fvcSmooth.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 int main(int argc, char *argv[])
56 {
57  #include "postProcess.H"
58 
59  #include "setRootCaseLists.H"
60  #include "createTime.H"
61  #include "createDynamicFvMesh.H"
62  #include "initContinuityErrs.H"
63  #include "createDyMControls.H"
64  #include "createFields.H"
65  #include "createUf.H"
66  #include "CourantNo.H"
67  #include "setInitialDeltaT.H"
68 
69  volScalarField& p = mixture.p();
70  volScalarField& T = mixture.T();
71  const volScalarField& psi1 = mixture.thermo1().psi();
72  const volScalarField& psi2 = mixture.thermo2().psi();
73 
74  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75  Info<< "\nStarting time loop\n" << endl;
76 
77  while (runTime.run())
78  {
79  #include "readDyMControls.H"
80 
81  // Store divU from the previous mesh so that it can be mapped
82  // and used in correctPhi to ensure the corrected phi has the
83  // same divergence
85 
86  if (LTS)
87  {
88  #include "setRDeltaT.H"
89  }
90  else
91  {
92  #include "CourantNo.H"
93  #include "alphaCourantNo.H"
94  #include "setDeltaT.H"
95  }
96 
97  runTime++;
98 
99  Info<< "Time = " << runTime.timeName() << nl << endl;
100 
101  // --- Pressure-velocity PIMPLE corrector loop
102  while (pimple.loop())
103  {
104  if (pimple.firstIter() || moveMeshOuterCorrectors)
105  {
106  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
107 
108  mesh.update();
109 
110  if (mesh.changing())
111  {
112 
113  MRF.update();
114 
115  Info<< "Execution time for mesh.update() = "
116  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
117  << " s" << endl;
118 
119  gh = (g & mesh.C()) - ghRef;
120  ghf = (g & mesh.Cf()) - ghRef;
121 
122  if (correctPhi)
123  {
124  // Calculate absolute flux
125  // from the mapped surface velocity
126  phi = mesh.Sf() & Uf;
127 
128  #include "correctPhi.H"
129 
130  // Make the fluxes relative to the mesh motion
132 
133  mixture.correct();
134  }
135 
136  if (checkMeshCourantNo)
137  {
138  #include "meshCourantNo.H"
139  }
140  }
141  }
142 
143  #include "alphaControls.H"
145 
146  turbulence.correctPhasePhi();
147 
148  #include "UEqn.H"
149  #include "TEqn.H"
150 
151  // --- Pressure corrector loop
152  while (pimple.correct())
153  {
154  #include "pEqn.H"
155  }
156 
157  if (pimple.turbCorr())
158  {
159  turbulence.correct();
160  }
161  }
162 
163  runTime.write();
164 
165  Info<< "ExecutionTime = "
166  << runTime.elapsedCpuTime()
167  << " s\n\n" << endl;
168  }
169 
170  Info<< "End\n" << endl;
171 
172  return 0;
173 }
174 
175 
176 // ************************************************************************* //
pimpleNoLoopControl & pimple
surfaceScalarField & phi
IOMRFZoneList & MRF
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:256
correctPhi
checkMeshCourantNo
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
Creates and initialises the velocity velocity field Uf.
autoPtr< surfaceVectorField > Uf
const surfaceScalarField & ghf
static const char nl
Definition: Ostream.H:265
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
moveMeshOuterCorrectors
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
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:188
U
Definition: pEqn.H:72
const volScalarField & gh
messageStream Info
zeroField divU
Definition: alphaSuSp.H:3
const dimensionedVector & g
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.