PDRFoam.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) 2011-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  PDRFoam
26 
27 Description
28  Solver for compressible premixed/partially-premixed combustion with
29  turbulence modelling.
30 
31  Combusting RANS code using the b-Xi two-equation model.
32  Xi may be obtained by either the solution of the Xi transport
33  equation or from an algebraic exression. Both approaches are
34  based on Gulder's flame speed correlation which has been shown
35  to be appropriate by comparison with the results from the
36  spectral model.
37 
38  Strain effects are incorporated directly into the Xi equation
39  but not in the algebraic approximation. Further work need to be
40  done on this issue, particularly regarding the enhanced removal rate
41  caused by flame compression. Analysis using results of the spectral
42  model will be required.
43 
44  For cases involving very lean Propane flames or other flames which are
45  very strain-sensitive, a transport equation for the laminar flame
46  speed is present. This equation is derived using heuristic arguments
47  involving the strain time scale and the strain-rate at extinction.
48  the transport velocity is the same as that for the Xi equation.
49 
50  For large flames e.g. explosions additional modelling for the flame
51  wrinkling due to surface instabilities may be applied.
52 
53  PDR (porosity/distributed resistance) modelling is included to handle
54  regions containing blockages which cannot be resolved by the mesh.
55 
56  The fields used by this solver are:
57 
58  betav: Volume porosity
59  Lobs: Average diameter of obstacle in cell (m)
60  Aw: Obstacle surface area per unit volume (1/m)
61  CR: Drag tensor (1/m)
62  CT: Turbulence generation parameter (1/m)
63  Nv: Number of obstacles in cell per unit volume (m^-2)
64  nsv: Tensor whose diagonal indicates the number to substract from
65  Nv to get the number of obstacles crossing the flow in each
66  direction.
67 
68 \*---------------------------------------------------------------------------*/
69 
70 #include "fvCFD.H"
71 #include "psiuReactionThermo.H"
73 #include "laminarFlameSpeed.H"
74 #include "XiModel.H"
75 #include "PDRDragModel.H"
76 #include "ignition.H"
77 #include "Switch.H"
78 #include "bound.H"
79 #include "pimpleControl.H"
80 
81 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83 int main(int argc, char *argv[])
84 {
85  #include "setRootCase.H"
86 
87  #include "createTime.H"
88  #include "createMesh.H"
89 
90  pimpleControl pimple(mesh);
91 
92  #include "readCombustionProperties.H"
93  #include "readGravitationalAcceleration.H"
94  #include "createFields.H"
95  #include "createMRF.H"
96  #include "initContinuityErrs.H"
97  #include "createTimeControls.H"
98  #include "compressibleCourantNo.H"
99  #include "setInitialDeltaT.H"
100 
101  scalar StCoNum = 0.0;
102 
103  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104 
105  Info<< "\nStarting time loop\n" << endl;
106 
107  while (runTime.run())
108  {
109  #include "createTimeControls.H"
110  #include "compressibleCourantNo.H"
111  #include "setDeltaT.H"
112 
113  runTime++;
114  Info<< "\n\nTime = " << runTime.timeName() << endl;
115 
116  #include "rhoEqn.H"
117 
118  // --- Pressure-velocity PIMPLE corrector loop
119  while (pimple.loop())
120  {
121  #include "UEqn.H"
122 
123  // --- Pressure corrector loop
124  while (pimple.correct())
125  {
126  #include "bEqn.H"
127  #include "ftEqn.H"
128  #include "EauEqn.H"
129  #include "EaEqn.H"
130 
131  if (!ign.ignited())
132  {
133  thermo.heu() == thermo.he();
134  }
135 
136  #include "pEqn.H"
137  }
138 
139  if (pimple.turbCorr())
140  {
141  turbulence->correct();
142  }
143  }
144 
145  runTime.write();
146 
147  Info<< "\nExecutionTime = "
148  << runTime.elapsedCpuTime()
149  << " s\n" << endl;
150  }
151 
152  Info<< "\n end\n";
153 
154  return 0;
155 }
156 
157 
158 // ************************************************************************* //
Bound the given scalar field if it has gone unbounded.
messageStream Info
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
dynamicFvMesh & mesh
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Read the control parameters used by setDeltaT.
int main(int argc, char *argv[])
Definition: postCalc.C:54
psiReactionThermo & thermo
Definition: createFields.H:32