overCompressibleInterDyMFoam.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) 2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  overCompressibleInterDyMFoam
28 
29 Group
30  grpMultiphaseSolvers
31 
32 Description
33  Solver for two compressible, non-isothermal, immiscible fluids using VOF
34  (i.e. volume of fluid) phase-fraction based interface capturing approach.
35 
36  This solver supports dynamic mesh motions including overset cases.
37 
38  The momentum and other fluid properties are of the "mixture" and a single
39  momentum equation is solved.
40 
41  Either mixture or two-phase transport modelling may be selected. In the
42  mixture approach, a single laminar, RAS or LES model is selected to model
43  the momentum stress. In the Euler-Euler two-phase approach separate
44  laminar, RAS or LES selected models are selected for each of the phases.
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #include "fvCFD.H"
49 #include "dynamicFvMesh.H"
50 #include "CMULES.H"
51 #include "EulerDdtScheme.H"
52 #include "CrankNicolsonDdtScheme.H"
53 #include "subCycle.H"
55 #include "pimpleControl.H"
56 #include "fvOptions.H"
57 #include "fvcSmooth.H"
58 
59 #include "cellCellStencilObject.H"
60 #include "localMin.H"
61 #include "interpolationCellPoint.H"
62 #include "transform.H"
63 #include "fvMeshSubset.H"
64 #include "oversetAdjustPhi.H"
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 int main(int argc, char *argv[])
69 {
70  argList::addNote
71  (
72  "Solver for two compressible, non-isothermal, immiscible fluids"
73  " using VOF phase-fraction based interface capturing approach.\n"
74  "Supports dynamic mesh motions including overset cases."
75  );
76 
77  #include "postProcess.H"
78 
79  #include "addCheckCaseOptions.H"
80  #include "setRootCaseLists.H"
81  #include "createTime.H"
82  #include "createDynamicFvMesh.H"
83  #include "createDyMControls.H"
84  #include "createFields.H"
85 
86  volScalarField& p = mixture.p();
87  volScalarField& T = mixture.T();
88  const volScalarField& psi1 = mixture.thermo1().psi();
89  const volScalarField& psi2 = mixture.thermo2().psi();
90 
91  #include "correctPhi.H"
92  #include "createUf.H"
93 
94  if (!LTS)
95  {
96  #include "CourantNo.H"
97  #include "setInitialDeltaT.H"
98  }
99 
100  #include "setCellMask.H"
101  #include "setInterpolatedCells.H"
102 
103  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104 
105  Info<< "\nStarting time loop\n" << endl;
106 
107  while (runTime.run())
108  {
109  #include "readDyMControls.H"
110 
111  if (LTS)
112  {
113  #include "setRDeltaT.H"
114  }
115  else
116  {
117  #include "CourantNo.H"
118  #include "alphaCourantNo.H"
119  #include "setDeltaT.H"
120  }
121 
122  ++runTime;
123 
124  Info<< "Time = " << runTime.timeName() << nl << endl;
125 
126  // --- Pressure-velocity PIMPLE corrector loop
127  while (pimple.loop())
128  {
129  if (pimple.firstIter() || moveMeshOuterCorrectors)
130  {
131  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
132 
133  mesh.update();
134 
135  if (mesh.changing())
136  {
137  Info<< "Execution time for mesh.update() = "
138  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
139  << " s" << endl;
140 
141  // Do not apply previous time-step mesh compression flux
142  // if the mesh topology changed
143  if (mesh.topoChanging())
144  {
145  talphaPhi1Corr0.clear();
146  }
147 
148  gh = (g & mesh.C()) - ghRef;
149  ghf = (g & mesh.Cf()) - ghRef;
150 
151  // Update cellMask field for blocking out hole cells
152  #include "setCellMask.H"
153  #include "setInterpolatedCells.H"
154  #include "correctPhiFaceMask.H"
155 
156  mixture.correct();
157 
158  // Make the flux relative to the mesh motion
160 
161  }
162 
163  if (mesh.changing() && checkMeshCourantNo)
164  {
165  #include "meshCourantNo.H"
166  }
167  }
168 
169  #include "alphaControls.H"
170  #include "compressibleAlphaEqnSubCycle.H"
171 
172  rhoPhi *= faceMask;
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 // ************************************************************************* //
Sets blocked cells mask field.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
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
checkMeshCourantNo
Sets blocked cells mask field.
const surfaceScalarField & ghf
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
dynamicFvMesh & mesh
const volScalarField & psi2
3D tensor transformation operations.
Creates and initialises the velocity velocity field Uf.
faceMask
Definition: setCellMask.H:43
const uniformDimensionedVectorField & g
bool LTS
Definition: createRDeltaT.H:1
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
moveMeshOuterCorrectors
const volScalarField & T
Correct phi on new faces C-I faces.
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
rhoPhi
Definition: rhoEqn.H:10
const volScalarField & gh
messageStream Info
Information stream (stdout output on master, null elsewhere)
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
Required Classes.
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Adjust the balance of fluxes on the faces between interpolated and calculated to obey continuity...
tmp< surfaceScalarField > talphaPhi1Corr0