MPPICInterFoam.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) 2016-2020 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  MPPICInterFoam
28 
29 Description
30  Solver for two incompressible, isothermal immiscible fluids using a VOF
31  (volume of fluid) phase-fraction based interface capturing approach.
32  The momentum and other fluid properties are of the "mixture" and a single
33  momentum equation is solved.
34 
35  It includes MRF and an MPPIC cloud.
36 
37  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "fvCFD.H"
42 #include "CMULES.H"
43 #include "EulerDdtScheme.H"
44 #include "localEulerDdtScheme.H"
45 #include "CrankNicolsonDdtScheme.H"
46 #include "subCycle.H"
47 
50 #include "pimpleControl.H"
51 #include "fvOptions.H"
52 #include "CorrectPhi.H"
53 #include "fvcSmooth.H"
54 
55 #include "basicKinematicCloud.H"
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 int main(int argc, char *argv[])
60 {
61  argList::addNote
62  (
63  "Solver for two incompressible, isothermal immiscible fluids using"
64  " VOF phase-fraction based interface capturing.\n"
65  "Includes MRF and an MPPIC cloud."
66  );
67 
68  #include "postProcess.H"
69 
70  #include "addCheckCaseOptions.H"
71  #include "setRootCaseLists.H"
72  #include "createTime.H"
73  #include "createMesh.H"
74  #include "createControl.H"
75  #include "createTimeControls.H"
76  #include "initContinuityErrs.H"
77  #include "createFields.H"
78  #include "createAlphaFluxes.H"
79  #include "createFvOptions.H"
80  #include "correctPhi.H"
81 
82  turbulence->validate();
83 
84  if (!LTS)
85  {
86  #include "readTimeControls.H"
87  #include "CourantNo.H"
88  #include "setInitialDeltaT.H"
89  }
90 
91  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
92 
93  Info<< "\nStarting time loop\n" << endl;
94 
95  while (runTime.run())
96  {
97  #include "readTimeControls.H"
98 
99  if (LTS)
100  {
101  #include "setRDeltaT.H"
102  }
103  else
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  Info<< "Evolving " << kinematicCloud.name() << endl;
115 
116  kinematicCloud.evolve();
117 
118  // Update continuous phase volume fraction field
119  alphac = max(1.0 - kinematicCloud.theta(), alphacMin);
120  alphac.correctBoundaryConditions();
121 
122  Info<< "Continuous phase-1 volume fraction = "
123  << alphac.weightedAverage(mesh.Vsc()).value()
124  << " Min(alphac) = " << min(alphac).value()
125  << " Max(alphac) = " << max(alphac).value()
126  << endl;
127 
128  alphacf = fvc::interpolate(alphac);
129  alphaRhoPhic = alphacf*rhoPhi;
130  alphaPhic = alphacf*phi;
131  alphacRho = alphac*rho;
132 
133  fvVectorMatrix cloudSU(kinematicCloud.SU(U));
134  volVectorField cloudVolSUSu
135  (
136  IOobject
137  (
138  "cloudVolSUSu",
139  runTime.timeName(),
140  mesh
141  ),
142  mesh,
143  dimensionedVector(cloudSU.dimensions()/dimVolume, Zero),
145  );
146 
147  cloudVolSUSu.primitiveFieldRef() = -cloudSU.source()/mesh.V();
148  cloudVolSUSu.correctBoundaryConditions();
149 
150  cloudSU.source() = vector::zero;
151 
152  // --- Pressure-velocity PIMPLE corrector loop
153  while (pimple.loop())
154  {
155  #include "alphaControls.H"
156  #include "alphaEqnSubCycle.H"
157 
158  mixture.correct();
159 
160  #include "UEqn.H"
161 
162  // --- Pressure corrector loop
163  while (pimple.correct())
164  {
165  #include "pEqn.H"
166  }
167 
168  if (pimple.turbCorr())
169  {
170  turbulence->correct();
171  }
172  }
173 
174  runTime.write();
175 
176  runTime.printExecutionTime(Info);
177  }
178 
179  Info<< "End\n" << endl;
180 
181  return 0;
182 }
183 
184 
185 // ************************************************************************* //
const word zeroGradientType
A zeroGradient patch field type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Read the control parameters used by setDeltaT.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
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
Required Classes.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
U
Definition: pEqn.H:72
rhoPhi
Definition: rhoEqn.H:10
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:40
messageStream Info
Information stream (stdout output on master, null elsewhere)
Execute application functionObjects to post-process existing results.
Required Classes.
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Read the control parameters used by setDeltaT.
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127