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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  DPMFoam
29 
30 Group
31  grpLagrangianSolvers
32 
33 Description
34  Transient solver for the coupled transport of a single kinematic particle
35  cloud including the effect of the volume fraction of particles on the
36  continuous phase.
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #include "fvCFD.H"
43 #include "pimpleControl.H"
44 
45 #ifdef MPPIC
46  #include "basicKinematicCloud.H"
47  #define basicKinematicTypeCloud basicKinematicCloud
48 #else
50  #define basicKinematicTypeCloud basicKinematicCollidingCloud
51 #endif
52 
53 int main(int argc, char *argv[])
54 {
55  argList::addNote
56  (
57  "Transient solver for the coupled transport of a"
58  " single kinematic particle cloud including the effect"
59  " of the volume fraction of particles on the continuous phase."
60  );
61  argList::addOption
62  (
63  "cloud",
64  "name",
65  "specify alternative cloud name. default is 'kinematicCloud'"
66  );
67 
68  #include "postProcess.H"
69 
70  #include "addCheckCaseOptions.H"
71  #include "setRootCaseLists.H"
72  #include "createTime.H"
73  #include "createMesh.H"
74  #include "createControl.H"
75  #include "createTimeControls.H"
76  #include "createFields.H"
77  #include "initContinuityErrs.H"
78 
79  Info<< "\nStarting time loop\n" << endl;
80 
81  while (runTime.run())
82  {
83  #include "readTimeControls.H"
84  #include "CourantNo.H"
85  #include "setDeltaT.H"
86 
87  ++runTime;
88 
89  Info<< "Time = " << runTime.timeName() << nl << endl;
90 
91  continuousPhaseTransport.correct();
92  muc = rhoc*continuousPhaseTransport.nu();
93 
94  Info<< "Evolving " << kinematicCloud.name() << endl;
95  kinematicCloud.evolve();
96 
97  // Update continuous phase volume fraction field
98  alphac = max(1.0 - kinematicCloud.theta(), alphacMin);
99  alphac.correctBoundaryConditions();
100  alphacf = fvc::interpolate(alphac);
101  alphaPhic = alphacf*phic;
102 
103  fvVectorMatrix cloudSU(kinematicCloud.SU(Uc));
104  volVectorField cloudVolSUSu
105  (
106  IOobject
107  (
108  "cloudVolSUSu",
109  runTime.timeName(),
110  mesh
111  ),
112  mesh,
113  dimensionedVector(cloudSU.dimensions()/dimVolume, Zero),
114  zeroGradientFvPatchVectorField::typeName
115  );
116 
117  cloudVolSUSu.primitiveFieldRef() = -cloudSU.source()/mesh.V();
118  cloudVolSUSu.correctBoundaryConditions();
119  cloudSU.source() = Zero;
120 
121 // cloudVolSUSu.primitiveFieldRef() =
122 // (cloudSU.diag()*Uc() - cloudSU.source())/mesh.V();
123 // cloudVolSUSu.correctBoundaryConditions();
124 // cloudSU.source() = cloudSU.diag()*Uc();
125 
126 
127  // --- Pressure-velocity PIMPLE corrector loop
128  while (pimple.loop())
129  {
130  #include "UcEqn.H"
131 
132  // --- PISO loop
133  while (pimple.correct())
134  {
135  #include "pEqn.H"
136  }
137 
138  if (pimple.turbCorr())
139  {
140  continuousPhaseTurbulence->correct();
141  }
142  }
143 
144  runTime.write();
145 
146  runTime.printExecutionTime(Info);
147  }
148 
149  Info<< "End\n" << endl;
150 
151  return 0;
152 }
153 
154 
155 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
singlePhaseTransportModel.H
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);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(dimless, Zero));const word kinematicCloudName(args.getOrDefault< word >"cloud", "kinematicCloud"));Info<< "Constructing kinematicCloud "<< kinematicCloudName<< endl;basicKinematicTypeCloud kinematicCloud(kinematicCloudName, rhoc, Uc, muc, g);scalar alphacMin(1.0 -(kinematicCloud.particleProperties().subDict("constantProperties") .get< scalar >"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< DPMIncompressibleTurbulenceModel< singlePhaseTransportModel > > continuousPhaseTurbulence(DPMIncompressibleTurbulenceModel< singlePhaseTransportModel >::New(alphac, Uc, alphaPhic, phic, continuousPhaseTransport))
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
basicKinematicCollidingCloud.H
DPMIncompressibleTurbulenceModel.H
pimpleControl.H
phic
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
addCheckCaseOptions.H
Required Classes.
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:47
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
basicKinematicCloud.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
postProcess.H
Execute application functionObjects to post-process existing results.
UcEqn.H
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
Foam::nl
constexpr char nl
Definition: Ostream.H:404
createTimeControls.H
Read the control parameters used by setDeltaT.
readTimeControls.H
Read the control parameters used by setDeltaT.
createMesh.H
Required Variables.
createTime.H
fvCFD.H
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.