compressibleInterDyMFoam.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-2015 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 Group
28  grpMultiphaseSolvers grpMovingMeshSolvers
29 
30 Description
31  Solver for 2 compressible, non-isothermal immiscible fluids using a VOF
32  (volume of fluid) phase-fraction based interface capturing approach,
33  with optional mesh motion and mesh topology changes including adaptive
34  re-meshing.
35 
36  The momentum and other fluid properties are of the "mixture" and a single
37  momentum equation is solved.
38 
39  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #include "fvCFD.H"
44 #include "dynamicFvMesh.H"
45 #include "MULES.H"
46 #include "subCycle.H"
47 #include "interfaceProperties.H"
48 #include "twoPhaseMixture.H"
49 #include "twoPhaseMixtureThermo.H"
51 #include "pimpleControl.H"
52 #include "CorrectPhi.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 int main(int argc, char *argv[])
58 {
59  #include "setRootCase.H"
60  #include "createTime.H"
61  #include "createDynamicFvMesh.H"
62  #include "initContinuityErrs.H"
63 
64  pimpleControl pimple(mesh);
65 
66  #include "createFields.H"
67  #include "createUf.H"
68  #include "createControls.H"
69  #include "CourantNo.H"
70  #include "setInitialDeltaT.H"
71 
72  turbulence->validate();
73 
74  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75  Info<< "\nStarting time loop\n" << endl;
76 
77  while (runTime.run())
78  {
79  #include "readControls.H"
80 
81  {
82  // Store divU from the previous mesh so that it can be mapped
83  // and used in correctPhi to ensure the corrected phi has the
84  // same divergence
86 
87  #include "CourantNo.H"
88  #include "setDeltaT.H"
89 
90  runTime++;
91 
92  Info<< "Time = " << runTime.timeName() << nl << endl;
93 
94  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
95 
96  // Do any mesh changes
97  mesh.update();
98 
99  if (mesh.changing())
100  {
101  Info<< "Execution time for mesh.update() = "
102  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
103  << " s" << endl;
104 
105  gh = (g & mesh.C()) - ghRef;
106  ghf = (g & mesh.Cf()) - ghRef;
107  }
108 
109  if (mesh.changing() && correctPhi)
110  {
111  // Calculate absolute flux from the mapped surface velocity
112  phi = mesh.Sf() & Uf;
113 
114  #include "correctPhi.H"
115 
116  // Make the fluxes relative to the mesh motion
118  }
119  }
120 
121  if (mesh.changing() && checkMeshCourantNo)
122  {
123  #include "meshCourantNo.H"
124  }
125 
126  turbulence->correct();
127 
128  // --- Pressure-velocity PIMPLE corrector loop
129  while (pimple.loop())
130  {
131  #include "alphaEqnsSubCycle.H"
132 
133  // correct interface on first PIMPLE corrector
134  if (pimple.corr() == 1)
135  {
136  interface.correct();
137  }
138 
140 
141  #include "UEqn.H"
142  #include "TEqn.H"
143 
144  // --- Pressure corrector loop
145  while (pimple.correct())
146  {
147  #include "pEqn.H"
148  }
149  }
150 
151  rho = alpha1*rho1 + alpha2*rho2;
152 
153  runTime.write();
154 
155  Info<< "ExecutionTime = "
156  << runTime.elapsedCpuTime()
157  << " s\n\n" << endl;
158  }
159 
160  Info<< "End\n" << endl;
161 
162  return 0;
163 }
164 
165 
166 // ************************************************************************* //
rhoPhi
rhoPhi
Definition: rhoEqn.H:10
ghf
const surfaceScalarField & ghf
Definition: setRegionFluidFields.H:35
subCycle.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
createUf.H
Creates and initialises the velocity velocity field Uf.
pimpleControl.H
U
U
Definition: pEqn.H:46
twoPhaseMixtureThermo.H
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:74
interface
interfaceProperties interface(alpha1, U, mixture())
solve
rhoEqn solve()
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
checkMeshCourantNo
bool checkMeshCourantNo
Definition: readControls.H:6
interfaceProperties.H
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
rho2
rho2
Definition: pEqn.H:125
alpha2
alpha2
Definition: alphaEqn.H:112
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
gh
const volScalarField & gh
Definition: setRegionFluidFields.H:34
rho1
rho1
Definition: pEqn.H:124
alpha1
volScalarField & alpha1
Definition: createFields.H:15
rho
rho
Definition: pEqn.H:3
setRootCase.H
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
CorrectPhi.H
createTime.H
dynamicFvMesh.H
fvCFD.H
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
fixedFluxPressureFvPatchScalarField.H
MULES.H
MULES: Multidimensional universal limiter for explicit solution.
Uf
Uf
Definition: pEqn.H:78
createDynamicFvMesh.H
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
turbulentFluidThermoModel.H
twoPhaseMixture.H
correctPhi
bool correctPhi
Definition: readControls.H:3