driftFluxFoam.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  driftFluxFoam
26 
27 Group
28  grpMultiphaseSolvers
29 
30 Description
31  Solver for 2 incompressible fluids using the mixture approach with the
32  drift-flux approximation for relative motion of the phases.
33 
34  Used for simulating the settling of the dispersed phase and other similar
35  separation problems.
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "fvCFD.H"
40 #include "CMULES.H"
41 #include "subCycle.H"
43 #include "relativeVelocityModel.H"
44 #include "turbulenceModel.H"
46 #include "pimpleControl.H"
47 #include "fvOptions.H"
49 #include "gaussLaplacianScheme.H"
50 #include "uncorrectedSnGrad.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 int main(int argc, char *argv[])
55 {
56  #include "setRootCase.H"
57 
58  #include "createTime.H"
59  #include "createMesh.H"
60 
61  pimpleControl pimple(mesh);
62 
63  #include "createTimeControls.H"
64  #include "createFields.H"
65  #include "createMRF.H"
66  #include "createFvOptions.H"
67  #include "initContinuityErrs.H"
68 
69  turbulence->validate();
70 
71  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 
73  Info<< "\nStarting time loop\n" << endl;
74 
75  while (runTime.run())
76  {
77  #include "readTimeControls.H"
78  #include "CourantNo.H"
79  #include "setDeltaT.H"
80 
81  runTime++;
82 
83  Info<< "Time = " << runTime.timeName() << nl << endl;
84 
85  // --- Pressure-velocity PIMPLE corrector loop
86  while (pimple.loop())
87  {
88  #include "alphaControls.H"
89 
90  UdmModel.correct();
91 
92  #include "alphaEqnSubCycle.H"
93 
94  mixture.correct();
95 
96  #include "UEqn.H"
97 
98  // --- Pressure corrector loop
99  while (pimple.correct())
100  {
101  #include "pEqn.H"
102  }
103 
104  if (pimple.turbCorr())
105  {
106  turbulence->correct();
107  }
108  }
109 
110  runTime.write();
111 
112  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
113  << " ClockTime = " << runTime.elapsedClockTime() << " s"
114  << nl << endl;
115  }
116 
117  Info<< "End\n" << endl;
118 
119  return 0;
120 }
121 
122 
123 // ************************************************************************* //
CMULES.H
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
fvOptions.H
subCycle.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
mixture
Info<< "Reading field p_rgh\n"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating phaseChangeTwoPhaseMixture\n"<< endl;autoPtr< phaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:33
pimpleControl.H
createFvOptions.H
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
incompressibleTwoPhaseInteractingMixture.H
relativeVelocityModel.H
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
CompressibleTurbulenceModel.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
setRootCase.H
createTimeControls.H
Read the control parameters used by setDeltaT.
readTimeControls.H
Read the control parameters used by setDeltaT.
createMesh.H
gaussLaplacianScheme.H
createTime.H
fvCFD.H
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
fixedFluxPressureFvPatchScalarField.H
UdmModel
relativeVelocityModel & UdmModel(UdmModelPtr())
turbulenceModel.H
uncorrectedSnGrad.H