DPMFoam.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) 2013-2018 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  DPMFoam
26 
27 Description
28  Transient solver for the coupled transport of a single kinematic particle
29  cloud including the effect of the volume fraction of particles on the
30  continuous phase.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #include "fvCFD.H"
37 #include "pimpleControl.H"
38 
39 #ifdef MPPIC
41  #define basicKinematicTypeCloud basicKinematicMPPICCloud
42 #else
44  #define basicKinematicTypeCloud basicKinematicCollidingCloud
45 #endif
46 
47 int main(int argc, char *argv[])
48 {
49  argList::addOption
50  (
51  "cloudName",
52  "name",
53  "specify alternative cloud name. default is 'kinematicCloud'"
54  );
55 
56  #include "postProcess.H"
57 
58  #include "setRootCaseLists.H"
59  #include "createTime.H"
60  #include "createMesh.H"
61  #include "createControl.H"
62  #include "createTimeControls.H"
63  #include "createFields.H"
64  #include "initContinuityErrs.H"
65 
66  Info<< "\nStarting time loop\n" << endl;
67 
68  while (runTime.run())
69  {
70  #include "readTimeControls.H"
71  #include "CourantNo.H"
72  #include "setDeltaT.H"
73 
74  runTime++;
75 
76  Info<< "Time = " << runTime.timeName() << nl << endl;
77 
78  continuousPhaseTransport.correct();
79  muc = rhoc*continuousPhaseTransport.nu();
80 
81  Info<< "Evolving " << kinematicCloud.name() << endl;
82  kinematicCloud.evolve();
83 
84  // Update continuous phase volume fraction field
85  alphac = max(1.0 - kinematicCloud.theta(), alphacMin);
86  alphac.correctBoundaryConditions();
87  alphacf = fvc::interpolate(alphac);
88  alphaPhic = alphacf*phic;
89 
90  fvVectorMatrix cloudSU(kinematicCloud.SU(Uc));
91  volVectorField cloudVolSUSu
92  (
93  IOobject
94  (
95  "cloudVolSUSu",
96  runTime.timeName(),
97  mesh
98  ),
99  mesh,
101  (
102  "0",
103  cloudSU.dimensions()/dimVolume,
104  Zero
105  ),
106  zeroGradientFvPatchVectorField::typeName
107  );
108 
109  cloudVolSUSu.primitiveFieldRef() = -cloudSU.source()/mesh.V();
110  cloudVolSUSu.correctBoundaryConditions();
111  cloudSU.source() = Zero;
112 
113  // --- Pressure-velocity PIMPLE corrector loop
114  while (pimple.loop())
115  {
116  #include "UcEqn.H"
117 
118  // --- PISO loop
119  while (pimple.correct())
120  {
121  #include "pEqn.H"
122  }
123 
124  if (pimple.turbCorr())
125  {
126  continuousPhaseTurbulence->correct();
127  }
128  }
129 
130  runTime.write();
131 
132  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
133  << " ClockTime = " << runTime.elapsedClockTime() << " s"
134  << nl << endl;
135  }
136 
137  Info<< "End\n" << endl;
138 
139  return 0;
140 }
141 
142 
143 // ************************************************************************* //
pimpleNoLoopControl & pimple
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Info<< "Reading field U\"<< endl;volVectorField Uc(IOobject(IOobject::groupName("U", continuousPhaseName), runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading/calculating continuous-phase face flux field phic\"<< endl;surfaceScalarField phic(IOobject(IOobject::groupName("phi", continuousPhaseName), runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), linearInterpolate(Uc) &mesh.Sf());label pRefCell=0;scalar pRefValue=0.0;setRefCell(p, pimple.dict(), pRefCell, pRefValue);mesh.setFluxRequired(p.name());Info<< "Creating turbulence model\"<< endl;singlePhaseTransportModel continuousPhaseTransport(Uc, phic);dimensionedScalar rhocValue(IOobject::groupName("rho", continuousPhaseName), dimDensity, continuousPhaseTransport.lookup(IOobject::groupName("rho", continuousPhaseName)));volScalarField rhoc(IOobject(rhocValue.name(), runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, rhocValue);volScalarField muc(IOobject(IOobject::groupName("mu", continuousPhaseName), runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), rhoc *continuousPhaseTransport.nu());Info<< "Creating field alphac\"<< endl;volScalarField alphac(IOobject(IOobject::groupName("alpha", continuousPhaseName), runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), mesh, dimensionedScalar("0", dimless, 0));word kinematicCloudName("kinematicCloud");args.optionReadIfPresent("cloudName", kinematicCloudName);Info<< "Constructing kinematicCloud "<< kinematicCloudName<< endl;basicKinematicTypeCloud kinematicCloud(kinematicCloudName, rhoc, Uc, muc, g);scalar alphacMin(1.0 - readScalar(kinematicCloud.particleProperties().subDict("constantProperties") .lookup("alphaMax")));alphac=max(1.0 - kinematicCloud.theta(), alphacMin);alphac.correctBoundaryConditions();surfaceScalarField alphacf("alphacf", fvc::interpolate(alphac));surfaceScalarField alphaPhic(IOobject::groupName("alphaPhi", continuousPhaseName), alphacf *phic);autoPtr< PhaseIncompressibleTurbulenceModel< singlePhaseTransportModel > > continuousPhaseTurbulence(PhaseIncompressibleTurbulenceModel< singlePhaseTransportModel >::New(alphac, Uc, alphaPhic, phic, continuousPhaseTransport))
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Read the control parameters used by setDeltaT.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
static const zero Zero
Definition: zero.H:97
static const char nl
Definition: Ostream.H:265
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
surfaceScalarField phic(mixture.cAlpha() *mag(phi/mesh.magSf()))
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
messageStream Info
Execute application functionObjects to post-process existing results.
Read the control parameters used by setDeltaT.