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  pimpleControl pimple(mesh);
84 
85  #include "createTimeControls.H"
86  #include "createDyMControls.H"
87  #include "createFields.H"
88 
89  volScalarField& p = mixture.p();
90  volScalarField& T = mixture.T();
91  const volScalarField& psi1 = mixture.thermo1().psi();
92  const volScalarField& psi2 = mixture.thermo2().psi();
93 
94  #include "correctPhi.H"
95  #include "createUf.H"
96 
97  if (!LTS)
98  {
99  #include "CourantNo.H"
100  #include "setInitialDeltaT.H"
101  }
102 
103  #include "setCellMask.H"
104  #include "setInterpolatedCells.H"
105 
106  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
107 
108  Info<< "\nStarting time loop\n" << endl;
109 
110  while (runTime.run())
111  {
112  #include "readControls.H"
113 
114  if (LTS)
115  {
116  #include "setRDeltaT.H"
117  }
118  else
119  {
120  #include "CourantNo.H"
121  #include "alphaCourantNo.H"
122  #include "setDeltaT.H"
123  }
124 
125  ++runTime;
126 
127  Info<< "Time = " << runTime.timeName() << nl << endl;
128 
129  // --- Pressure-velocity PIMPLE corrector loop
130  while (pimple.loop())
131  {
132  if (pimple.firstIter() || moveMeshOuterCorrectors)
133  {
134  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
135 
136  mesh.update();
137 
138  if (mesh.changing())
139  {
140  Info<< "Execution time for mesh.update() = "
141  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
142  << " s" << endl;
143 
144  // Do not apply previous time-step mesh compression flux
145  // if the mesh topology changed
146  if (mesh.topoChanging())
147  {
148  talphaPhi1Corr0.clear();
149  }
150 
151  gh = (g & mesh.C()) - ghRef;
152  ghf = (g & mesh.Cf()) - ghRef;
153 
154  // Update cellMask field for blocking out hole cells
155  #include "setCellMask.H"
156  #include "setInterpolatedCells.H"
157 
158  faceMask =
159  localMin<scalar>(mesh).interpolate(cellMask.oldTime());
160 
161  // Zero Uf on old faceMask (H-I)
162  Uf *= faceMask;
163 
165  // Update Uf and phi on new C-I faces
166  Uf += (1-faceMask)*Uint;
167 
168  // Update Uf boundary
169  forAll(Uf.boundaryField(), patchI)
170  {
171  Uf.boundaryFieldRef()[patchI] =
172  Uint.boundaryField()[patchI];
173  }
174 
175  phi = mesh.Sf() & Uf;
176 
177  // Correct phi on individual regions
178  if (correctPhi)
179  {
180  #include "correctPhi.H"
181  }
182 
183  mixture.correct();
184 
185  // Zero phi on current H-I
186  faceMask = localMin<scalar>(mesh).interpolate(cellMask);
187 
188  phi *= faceMask;
189  U *= cellMask;
190 
191  // Make the flux relative to the mesh motion
193 
194  }
195 
196  if (mesh.changing() && checkMeshCourantNo)
197  {
198  #include "meshCourantNo.H"
199  }
200  }
201 
202  #include "alphaControls.H"
203  #include "compressibleAlphaEqnSubCycle.H"
204 
206  (
207  localMin<scalar>(mesh).interpolate(cellMask)
208  );
209  rhoPhi *= faceMask;
210 
211  turbulence.correctPhasePhi();
212 
213  #include "UEqn.H"
214  volScalarField divUp("divUp", fvc::div(fvc::absolute(phi, U), p));
215  #include "TEqn.H"
216 
217  // --- Pressure corrector loop
218  while (pimple.correct())
219  {
220  #include "pEqn.H"
221  }
222 
223  if (pimple.turbCorr())
224  {
225  turbulence.correct();
226  }
227  }
228 
229  runTime.write();
230 
231  runTime.printExecutionTime(Info);
232  }
233 
234  Info<< "End\n" << endl;
235 
236  return 0;
237 }
238 
239 
240 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
rhoPhi
rhoPhi
Definition: rhoEqn.H:10
p
volScalarField & p
Definition: createFieldRefs.H:8
CMULES.H
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
ghf
const surfaceScalarField & ghf
Definition: setRegionFluidFields.H:18
fvOptions.H
turbulence
Info<< "Reading field U\n"<< 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\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
fvMeshSubset.H
subCycle.H
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:27
faceMask
surfaceScalarField faceMask(localMin< scalar >(mesh).interpolate(cellMask))
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
createUf.H
Creates and initialises the velocity velocity field Uf.
setCellMask.H
Sets blocked cells mask field.
setInterpolatedCells.H
Sets blocked cells mask field.
correctPhi
correctPhi
Definition: readDyMControls.H:3
gh
const volScalarField & gh
Definition: setRegionFluidFields.H:17
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
psi2
const volScalarField & psi2
Definition: setRegionFluidFields.H:31
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
interpolationCellPoint.H
pimpleControl.H
compressibleInterPhaseTransportModel.H
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Definition: fvcMeshPhi.C:70
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
addCheckCaseOptions.H
Required Classes.
Foam::Info
messageStream Info
EulerDdtScheme.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
CrankNicolsonDdtScheme.H
LTS
bool LTS
Definition: createRDeltaT.H:1
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
postProcess.H
Execute application functionObjects to post-process existing results.
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
psi1
const volScalarField & psi1
Definition: setRegionFluidFields.H:28
T
const volScalarField & T
Definition: createFieldRefs.H:2
U
U
Definition: pEqn.H:72
localMin.H
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:50
createTimeControls.H
Read the control parameters used by setDeltaT.
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
createTime.H
dynamicFvMesh.H
fvCFD.H
checkMeshCourantNo
checkMeshCourantNo
Definition: readDyMControls.H:9
transform.H
3D tensor transformation operations.
Foam::surfaceVectorField
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Definition: surfaceFieldsFwd.H:55
moveMeshOuterCorrectors
moveMeshOuterCorrectors
Definition: readDyMControls.H:15
cellCellStencilObject.H
mixture
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
talphaPhi1Corr0
tmp< surfaceScalarField > talphaPhi1Corr0
Definition: createAlphaFluxes.H:26
oversetAdjustPhi.H
Adjust the balance of fluxes on the faces between interpolated and calculated to obey continuity.
createDynamicFvMesh.H
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Definition: fvcMeshPhi.C:183
fvcSmooth.H
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...