DPMFoam.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) 2013-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  DPMFoam
26 
27 Group
28  grpLagrangianSolvers
29 
30 Description
31  Transient solver for the coupled transport of a single kinematic particle
32  cloud including the effect of the volume fraction of particles on the
33  continuous phase.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "fvCFD.H"
40 #include "pimpleControl.H"
42 
43 #ifdef MPPIC
45  #define basicKinematicTypeCloud basicKinematicMPPICCloud
46 #else
48  #define basicKinematicTypeCloud basicKinematicCollidingCloud
49 #endif
50 
51 int main(int argc, char *argv[])
52 {
53  argList::addOption
54  (
55  "cloudName",
56  "name",
57  "specify alternative cloud name. default is 'kinematicCloud'"
58  );
59 
60  #include "setRootCase.H"
61  #include "createTime.H"
62  #include "createMesh.H"
63 
64  pimpleControl pimple(mesh);
65 
66  #include "createTimeControls.H"
67  #include "readGravitationalAcceleration.H"
68  #include "createFields.H"
69  #include "initContinuityErrs.H"
70 
71  Info<< "\nStarting time loop\n" << endl;
72 
73  while (runTime.run())
74  {
75  #include "readTimeControls.H"
76  #include "CourantNo.H"
77  #include "setDeltaT.H"
78 
79  runTime++;
80 
81  Info<< "Time = " << runTime.timeName() << nl << endl;
82 
83  continuousPhaseTransport.correct();
84  muc = rhoc*continuousPhaseTransport.nu();
85 
86  Info<< "Evolving " << kinematicCloud.name() << endl;
87  kinematicCloud.evolve();
88 
89  // Update continuous phase volume fraction field
90  alphac = max(1.0 - kinematicCloud.theta(), alphacMin);
91  alphac.correctBoundaryConditions();
92  alphacf = fvc::interpolate(alphac);
93  alphaPhic = alphacf*phic;
94 
95  fvVectorMatrix cloudSU(kinematicCloud.SU(Uc));
96  volVectorField cloudVolSUSu
97  (
98  IOobject
99  (
100  "cloudVolSUSu",
101  runTime.timeName(),
102  mesh
103  ),
104  mesh,
106  (
107  "0",
108  cloudSU.dimensions()/dimVolume,
109  vector::zero
110  ),
111  zeroGradientFvPatchVectorField::typeName
112  );
113 
114  cloudVolSUSu.internalField() = -cloudSU.source()/mesh.V();
115  cloudVolSUSu.correctBoundaryConditions();
116  cloudSU.source() = vector::zero;
117 
118  // --- Pressure-velocity PIMPLE corrector loop
119  while (pimple.loop())
120  {
121  #include "UcEqn.H"
122 
123  // --- PISO loop
124  while (pimple.correct())
125  {
126  #include "pEqn.H"
127  }
128 
129  if (pimple.turbCorr())
130  {
131  continuousPhaseTurbulence->correct();
132  }
133  }
134 
135  runTime.write();
136 
137  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
138  << " ClockTime = " << runTime.elapsedClockTime() << " s"
139  << nl << endl;
140  }
141 
142  Info<< "End\n" << endl;
143 
144  return 0;
145 }
146 
147 
148 // ************************************************************************* //
basicKinematicMPPICCloud.H
singlePhaseTransportModel.H
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
basicKinematicCollidingCloud.H
PhaseIncompressibleTurbulenceModel.H
pimpleControl.H
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:48
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
continuousPhaseTurbulence
Info<< "Reading field U\n"<< endl;volVectorField Uc(IOobject(IOobject::groupName("U", continuousPhaseName), runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field p\n"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading/calculating continuous-phase face flux field phic\n"<< 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\n"<< 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\n"<< 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("alphaPhic", alphacf *phic);autoPtr< PhaseIncompressibleTurbulenceModel< singlePhaseTransportModel > > continuousPhaseTurbulence(PhaseIncompressibleTurbulenceModel< singlePhaseTransportModel >::New(alphac, Uc, alphaPhic, phic, continuousPhaseTransport))
Definition: createFields.H:156
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
UcEqn.H
phic
phic
Definition: alphaEqnsSubCycle.H:5
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
setRootCase.H
createTimeControls.H
Read the control parameters used by setDeltaT.
readTimeControls.H
Read the control parameters used by setDeltaT.
createMesh.H
createTime.H
fvCFD.H
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
fixedFluxPressureFvPatchScalarField.H
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58