DPMDyMFoam.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) 2017 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  DPMDyMFoam
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, with optional mesh motion and mesh topology changes.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "fvCFD.H"
38 #include "dynamicFvMesh.H"
41 #include "pimpleControl.H"
42 #include "CorrectPhi.H"
43 
44 #ifdef MPPIC
45  #include "basicKinematicCloud.H"
46  #define basicKinematicTypeCloud basicKinematicCloud
47 #else
49  #define basicKinematicTypeCloud basicKinematicCollidingCloud
50 #endif
51 
52 int main(int argc, char *argv[])
53 {
54  argList::addNote
55  (
56  "Transient solver for the coupled transport of a"
57  " single kinematic particle cloud including the effect"
58  " of the volume fraction of particles on the continuous phase.\n"
59  "With optional mesh motion and mesh topology changes."
60  );
61  argList::addOption
62  (
63  "cloudName",
64  "name",
65  "specify alternative cloud name. default is 'kinematicCloud'"
66  );
67 
68  #include "postProcess.H"
69 
70  #include "setRootCaseLists.H"
71  #include "createTime.H"
72  #include "createDynamicFvMesh.H"
73  #include "createDyMControls.H"
74  #include "createFields.H"
75  #include "createUcf.H"
76  #include "initContinuityErrs.H"
77 
78  Info<< "\nStarting time loop\n" << endl;
79 
80  while (runTime.run())
81  {
82  #include "readDyMControls.H"
83  #include "CourantNo.H"
84  #include "setDeltaT.H"
85 
86  ++runTime;
87 
88  Info<< "Time = " << runTime.timeName() << nl << endl;
89 
90  // Store the particle positions
91  kinematicCloud.storeGlobalPositions();
92 
93  mesh.update();
94 
95  // Calculate absolute flux from the mapped surface velocity
96  phic = mesh.Sf() & Ucf;
97 
98  if (mesh.changing() && correctPhi)
99  {
100  #include "correctPhic.H"
101  }
102 
103  // Make the flux relative to the mesh motion
104  fvc::makeRelative(phic, Uc);
105 
106  if (mesh.changing() && checkMeshCourantNo)
107  {
108  #include "meshCourantNo.H"
109  }
110 
111  continuousPhaseTransport.correct();
112  muc = rhoc*continuousPhaseTransport.nu();
113 
114  kinematicCloud.evolve();
115 
116  // Update continuous phase volume fraction field
117  alphac = max(1.0 - kinematicCloud.theta(), alphacMin);
118  alphac.correctBoundaryConditions();
119  alphacf = fvc::interpolate(alphac);
120  alphaPhic = alphacf*phic;
121 
122  fvVectorMatrix cloudSU(kinematicCloud.SU(Uc));
123  volVectorField cloudVolSUSu
124  (
125  IOobject
126  (
127  "cloudVolSUSu",
128  runTime.timeName(),
129  mesh
130  ),
131  mesh,
132  dimensionedVector(cloudSU.dimensions()/dimVolume, Zero),
134  );
135 
136  cloudVolSUSu.primitiveFieldRef() = -cloudSU.source()/mesh.V();
137  cloudVolSUSu.correctBoundaryConditions();
138  cloudSU.source() = Zero;
139 
140  // --- Pressure-velocity PIMPLE corrector loop
141  while (pimple.loop())
142  {
143  #include "UcEqn.H"
144 
145  // --- PISO loop
146  while (pimple.correct())
147  {
148  #include "pEqn.H"
149  }
150 
151  if (pimple.turbCorr())
152  {
153  continuousPhaseTurbulence->correct();
154  }
155  }
156 
157  runTime.write();
158 
159  runTime.printExecutionTime(Info);
160  }
161 
162  Info<< "End\n" << endl;
163 
164  return 0;
165 }
166 
167 
168 // ************************************************************************* //
const word zeroGradientType
A zeroGradient patch field type.
Creates and initialises the velocity velocity field Ucf.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
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.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
correctPhi
checkMeshCourantNo
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
pimpleControl & pimple
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
Calculates and outputs the mean and maximum Courant Numbers.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:40
messageStream Info
Information stream (stdout output on master, null elsewhere)
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);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, 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))
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:70
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133