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, with optional mesh motion and mesh topology changes.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #include "fvCFD.H"
35 #include "dynamicFvMesh.H"
38 #include "pimpleControl.H"
39 #include "CorrectPhi.H"
40 
41 #ifdef MPPIC
43  #define basicKinematicTypeCloud basicKinematicMPPICCloud
44 #else
46  #define basicKinematicTypeCloud basicKinematicCollidingCloud
47 #endif
48 
49 int main(int argc, char *argv[])
50 {
51  argList::addOption
52  (
53  "cloudName",
54  "name",
55  "specify alternative cloud name. default is 'kinematicCloud'"
56  );
57 
58  #include "postProcess.H"
59 
60  #include "setRootCaseLists.H"
61  #include "createTime.H"
62  #include "createDynamicFvMesh.H"
63  #include "createDyMControls.H"
64  #include "createFields.H"
65  #include "createUcfIfPresent.H"
66  #include "initContinuityErrs.H"
67 
68  Info<< "\nStarting time loop\n" << endl;
69 
70  while (runTime.run())
71  {
72  #include "readDyMControls.H"
73  #include "CourantNo.H"
74  #include "setDeltaT.H"
75 
76  runTime++;
77 
78  Info<< "Time = " << runTime.timeName() << nl << endl;
79 
80  // Store the particle positions
81  kinematicCloud.storeGlobalPositions();
82 
83  mesh.update();
84 
85  if (mesh.changing())
86  {
87  if (correctPhi)
88  {
89  // Calculate absolute flux from the mapped surface velocity
90  phic = mesh.Sf() & Ucf();
91 
92  #include "correctPhic.H"
93 
94  // Make the flux relative to the mesh motion
95  fvc::makeRelative(phic, Uc);
96  }
97 
99  {
100  #include "meshCourantNo.H"
101  }
102  }
103 
104  continuousPhaseTransport.correct();
105  muc = rhoc*continuousPhaseTransport.nu();
106 
107  Info<< "Evolving " << kinematicCloud.name() << endl;
108  kinematicCloud.evolve();
109 
110  // Update continuous phase volume fraction field
111  alphac = max(1.0 - kinematicCloud.theta(), alphacMin);
112  alphac.correctBoundaryConditions();
113  alphacf = fvc::interpolate(alphac);
114  alphaPhic = alphacf*phic;
115 
116  fvVectorMatrix cloudSU(kinematicCloud.SU(Uc));
117  volVectorField cloudVolSUSu
118  (
119  IOobject
120  (
121  "cloudVolSUSu",
122  runTime.timeName(),
123  mesh
124  ),
125  mesh,
127  (
128  "0",
129  cloudSU.dimensions()/dimVolume,
130  Zero
131  ),
132  zeroGradientFvPatchVectorField::typeName
133  );
134 
135  cloudVolSUSu.primitiveFieldRef() = -cloudSU.source()/mesh.V();
136  cloudVolSUSu.correctBoundaryConditions();
137  cloudSU.source() = Zero;
138 
139  // --- Pressure-velocity PIMPLE corrector loop
140  while (pimple.loop())
141  {
142  #include "UcEqn.H"
143 
144  // --- PISO loop
145  while (pimple.correct())
146  {
147  #include "pEqn.H"
148  }
149 
150  if (pimple.turbCorr())
151  {
152  continuousPhaseTurbulence->correct();
153  }
154  }
155 
156  runTime.write();
157 
158  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
159  << " ClockTime = " << runTime.elapsedClockTime() << " s"
160  << nl << endl;
161  }
162 
163  Info<< "End\n" << endl;
164 
165  return 0;
166 }
167 
168 
169 // ************************************************************************* //
pimpleNoLoopControl & pimple
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
autoPtr< surfaceVectorField > Ucf
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
correctPhi
checkMeshCourantNo
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(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))
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Creates and initialises the continuous phase face velocity field Ufc if required. ...
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.
Calculates and outputs the mean and maximum Courant Numbers.
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.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75