driftFluxFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 "momentumTransportModel.H"
42 #include "pimpleControl.H"
43 #include "pressureReference.H"
44 #include "fvModels.H"
45 #include "fvConstraints.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 int main(int argc, char *argv[])
50 {
51  #include "postProcess.H"
52 
53  #include "setRootCaseLists.H"
54  #include "createTime.H"
55  #include "createMesh.H"
56  #include "createControl.H"
57  #include "createTimeControls.H"
58  #include "createFields.H"
59  #include "initContinuityErrs.H"
60 
61  volScalarField& alpha2(mixture.alpha2());
62  const dimensionedScalar& rho1 = mixture.rhod();
63  const dimensionedScalar& rho2 = mixture.rhoc();
64 
65  turbulence->validate();
66 
67  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69  Info<< "\nStarting time loop\n" << endl;
70 
71  while (pimple.run(runTime))
72  {
73  #include "readTimeControls.H"
74  #include "CourantNo.H"
75  #include "setDeltaT.H"
76 
77  runTime++;
78 
79  Info<< "Time = " << runTime.userTimeName() << nl << endl;
80 
81  // --- Pressure-velocity PIMPLE corrector loop
82  while (pimple.loop())
83  {
84  mixture.correct();
85 
86  fvModels.correct();
87 
88  #include "alphaEqnSubCycle.H"
89 
90  mixture.correct();
91 
92  #include "UEqn.H"
93 
94  // --- Pressure corrector loop
95  while (pimple.correct())
96  {
97  #include "pEqn.H"
98  }
99 
100  if (pimple.turbCorr())
101  {
102  turbulence->correct();
103  }
104  }
105 
106  runTime.write();
107 
108  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
109  << " ClockTime = " << runTime.elapsedClockTime() << " s"
110  << nl << endl;
111  }
112 
113  Info<< "End\n" << endl;
114 
115  return 0;
116 }
117 
118 
119 // ************************************************************************* //
pimpleNoLoopControl & pimple
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:353
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Read the control parameters used by setDeltaT.
alpha2
Definition: alphaEqn.H:115
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
static const char nl
Definition: Ostream.H:260
Foam::fvModels & fvModels
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.