reactingParcelFoam.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) 2011-2020 OpenFOAM Foundation
9 Copyright (C) 2018-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Application
28 reactingParcelFoam
29
30Group
31 grpLagrangianSolvers
32
33Description
34 Transient solver for compressible, turbulent flow with a reacting,
35 multiphase particle cloud, and surface film modelling.
36
37\*---------------------------------------------------------------------------*/
38
39#include "fvCFD.H"
40#include "dynamicFvMesh.H"
42#include "surfaceFilmModel.H"
43#include "rhoReactionThermo.H"
44#include "CombustionModel.H"
45#include "radiationModel.H"
46#include "SLGThermo.H"
47#include "fvOptions.H"
48#include "pimpleControl.H"
49#include "pressureControl.H"
50#include "CorrectPhi.H"
51#include "localEulerDdtScheme.H"
52#include "fvcSmooth.H"
53#include "cloudMacros.H"
54
55#ifndef CLOUD_BASE_TYPE
56 #define CLOUD_BASE_TYPE ReactingMultiphase
57 #define CLOUD_BASE_TYPE_NAME "reacting"
58#endif
59
60#include CLOUD_INCLUDE_FILE(CLOUD_BASE_TYPE)
61#define basicReactingTypeCloud CLOUD_TYPE(CLOUD_BASE_TYPE)
62
63// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64
65int main(int argc, char *argv[])
66{
67 argList::addNote
68 (
69 "Transient solver for compressible, turbulent flow"
70 " with reacting, multiphase particle clouds"
71 " and surface film modelling."
72 );
73
74 #define CREATE_MESH createMeshesPostProcess.H
75 #include "postProcess.H"
76
77 #include "addCheckCaseOptions.H"
78 #include "setRootCaseLists.H"
79 #include "createTime.H"
80 #include "createDynamicFvMesh.H"
81 #include "createDyMControls.H"
82 #include "createFields.H"
83 #include "createFieldRefs.H"
84 #include "createRegionControls.H"
85 #include "initContinuityErrs.H"
86 #include "createRhoUfIfPresent.H"
87
88 turbulence->validate();
89
90 if (!LTS)
91 {
92 #include "compressibleCourantNo.H"
93 #include "setInitialDeltaT.H"
94 }
95
96 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97
98 Info<< "\nStarting time loop\n" << endl;
99
100 while (runTime.run())
101 {
102 #include "readDyMControls.H"
103
104 // Store divrhoU from the previous mesh
105 // so that it can be mapped and used in correctPhi
106 // to ensure the corrected phi has the same divergence
107 autoPtr<volScalarField> divrhoU;
109 {
110 divrhoU.reset
111 (
112 new volScalarField
113 (
114 "divrhoU",
115 fvc::div(fvc::absolute(phi, rho, U))
116 )
117 );
118 }
119
120 if (LTS)
121 {
122 #include "setRDeltaT.H"
123 }
124 else
125 {
126 #include "compressibleCourantNo.H"
127 #include "setMultiRegionDeltaT.H"
128 }
129
130 ++runTime;
131
132 Info<< "Time = " << runTime.timeName() << nl << endl;
133
134 // Store momentum to set rhoUf for introduced faces.
135 autoPtr<volVectorField> rhoU;
136 if (solvePrimaryRegion && rhoUf.valid())
137 {
138 rhoU.reset(new volVectorField("rhoU", rho*U));
139 }
140
141 // Store the particle positions
142 parcels.storeGlobalPositions();
143
144 // Do any mesh changes
145 mesh.update();
146
147 if (solvePrimaryRegion && mesh.changing())
148 {
149 gh = (g & mesh.C()) - ghRef;
150 ghf = (g & mesh.Cf()) - ghRef;
151
152 MRF.update();
153
154 if (correctPhi)
155 {
156 // Calculate absolute flux
157 // from the mapped surface velocity
158 phi = mesh.Sf() & rhoUf();
159
160 #include "../../compressible/rhoPimpleFoam/correctPhi.H"
161
162 // Make the fluxes relative to the mesh-motion
163 fvc::makeRelative(phi, rho, U);
164 }
165
167 {
168 #include "meshCourantNo.H"
169 }
170 }
171
172 parcels.evolve();
173 surfaceFilm.evolve();
174
176 {
177 if (pimple.nCorrPIMPLE() <= 1)
178 {
179 #include "rhoEqn.H"
180 }
181
182 // --- PIMPLE loop
183 while (pimple.loop())
184 {
185 #include "UEqn.H"
186 #include "YEqn.H"
187 #include "EEqn.H"
188
189 // --- Pressure corrector loop
190 while (pimple.correct())
191 {
192 #include "pEqn.H"
193 }
194
195 if (pimple.turbCorr())
196 {
197 turbulence->correct();
198 }
199 }
200
201 rho = thermo.rho();
202 }
203
204 runTime.write();
205
206 runTime.printExecutionTime(Info);
207 }
208
209 Info<< "End\n" << endl;
210
211 return 0;
212}
213
214
215// ************************************************************************* //
bool solvePrimaryRegion(pimpleDict.getOrDefault("solvePrimaryRegion", true))
Required Classes.
const uniformDimensionedVectorField & g
surfaceScalarField & phi
const surfaceScalarField & ghf
IOMRFZoneList & MRF
const volScalarField & gh
pimpleControl & pimple
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
C-preprocessor cloud macros.
U
Definition: pEqn.H:72
regionModels::surfaceFilmModel & surfaceFilm
dynamicFvMesh & mesh
engineTime & runTime
bool LTS
Definition: createRDeltaT.H:1
Creates and initialises the velocity field rhoUf if required.
autoPtr< surfaceVectorField > rhoUf
compressible::turbulenceModel & turbulence
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Calculates and outputs the mean and maximum Courant Numbers.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
Execute application functionObjects to post-process existing results.
checkMeshCourantNo
correctPhi