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-2016 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 Description
28  Solver for 2 incompressible fluids using the mixture approach with the
29  drift-flux approximation for relative motion of the phases.
30 
31  Used for simulating the settling of the dispersed phase and other similar
32  separation problems.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "fvCFD.H"
37 #include "CMULES.H"
38 #include "subCycle.H"
40 #include "relativeVelocityModel.H"
41 #include "turbulenceModel.H"
43 #include "pimpleControl.H"
44 #include "fvOptions.H"
45 #include "gaussLaplacianScheme.H"
46 #include "uncorrectedSnGrad.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 int main(int argc, char *argv[])
51 {
52  #include "postProcess.H"
53 
54  #include "setRootCase.H"
55  #include "createTime.H"
56  #include "createMesh.H"
57  #include "createControl.H"
58  #include "createTimeControls.H"
59  #include "createFields.H"
60  #include "createFvOptions.H"
61  #include "initContinuityErrs.H"
62 
63  volScalarField& alpha2(mixture.alpha2());
64  const dimensionedScalar& rho1 = mixture.rhod();
65  const dimensionedScalar& rho2 = mixture.rhoc();
66  relativeVelocityModel& UdmModel(UdmModelPtr());
67 
68  turbulence->validate();
69 
70  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72  Info<< "\nStarting time loop\n" << endl;
73 
74  while (runTime.run())
75  {
76  #include "readTimeControls.H"
77  #include "CourantNo.H"
78  #include "setDeltaT.H"
79 
80  runTime++;
81 
82  Info<< "Time = " << runTime.timeName() << nl << endl;
83 
84  // --- Pressure-velocity PIMPLE corrector loop
85  while (pimple.loop())
86  {
87  #include "alphaControls.H"
88 
89  UdmModel.correct();
90 
91  #include "alphaEqnSubCycle.H"
92 
93  mixture.correct();
94 
95  #include "UEqn.H"
96 
97  // --- Pressure corrector loop
98  while (pimple.correct())
99  {
100  #include "pEqn.H"
101  }
102 
103  if (pimple.turbCorr())
104  {
105  turbulence->correct();
106  }
107  }
108 
109  runTime.write();
110 
111  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
112  << " ClockTime = " << runTime.elapsedClockTime() << " s"
113  << nl << endl;
114  }
115 
116  Info<< "End\n" << endl;
117 
118  return 0;
119 }
120 
121 
122 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
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<< "Reading incompressibleTwoPhaseInteractingMixture\n"<< endl;incompressibleTwoPhaseInteractingMixture mixture(U, phi);volScalarField &alpha1(mixture.alpha1());volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mixture.rho());surfaceScalarField rhoPhi(IOobject("rhoPhi", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), fvc::interpolate(rho)*phi);autoPtr< relativeVelocityModel > UdmModelPtr(relativeVelocityModel::New(mixture, mixture))
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Read the control parameters used by setDeltaT.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
alpha2
Definition: alphaEqn.H:112
static const char nl
Definition: Ostream.H:262
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
Execute application functionObjects to post-process existing results.
Read the control parameters used by setDeltaT.
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.