PDRFoam.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) 2011-2021 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 expression. 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 subtract 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"
74 #include "laminarFlameSpeed.H"
75 #include "XiModel.H"
76 #include "PDRDragModel.H"
77 #include "ignition.H"
78 #include "Switch.H"
79 #include "bound.H"
80 #include "pimpleControl.H"
81 #include "fvModels.H"
82 #include "fvConstraints.H"
83 
84 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85 
86 int main(int argc, char *argv[])
87 {
88  #include "postProcess.H"
89 
90  #include "setRootCaseLists.H"
91  #include "createTime.H"
92  #include "createMesh.H"
93  #include "createControl.H"
94  #include "readCombustionProperties.H"
95  #include "readGravitationalAcceleration.H"
96  #include "createFields.H"
97  #include "createFieldRefs.H"
98  #include "initContinuityErrs.H"
99  #include "createTimeControls.H"
100  #include "compressibleCourantNo.H"
101  #include "setInitialDeltaT.H"
102 
103  turbulence->validate();
104  scalar StCoNum = 0.0;
105 
106  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
107 
108  Info<< "\nStarting time loop\n" << endl;
109 
110  while (pimple.run(runTime))
111  {
112  #include "readTimeControls.H"
113  #include "compressibleCourantNo.H"
114  #include "setDeltaT.H"
115 
116  runTime++;
117  Info<< "\n\nTime = " << runTime.timeName() << endl;
118 
119  #include "rhoEqn.H"
120 
121  // --- Pressure-velocity PIMPLE corrector loop
122  while (pimple.loop())
123  {
124  fvModels.correct();
125 
126  #include "UEqn.H"
127 
128  // --- Pressure corrector loop
129  while (pimple.correct())
130  {
131  #include "bEqn.H"
132  #include "ftEqn.H"
133  #include "EauEqn.H"
134  #include "EaEqn.H"
135 
136  if (!ign.ignited())
137  {
138  thermo.heu() == thermo.he();
139  }
140 
141  #include "pEqn.H"
142  }
143 
144  if (pimple.turbCorr())
145  {
146  turbulence->correct();
147  thermophysicalTransport->correct();
148  }
149  }
150 
151  runTime.write();
152 
153  Info<< "\nExecutionTime = "
154  << runTime.elapsedCpuTime()
155  << " s\n" << endl;
156  }
157 
158  Info<< "\n end\n";
159 
160  return 0;
161 }
162 
163 
164 // ************************************************************************* //
fluidReactionThermo & thermo
Definition: createFields.H:28
pimpleNoLoopControl & pimple
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:316
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Read the control parameters used by setDeltaT.
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Bound the given scalar field if it has gone unbounded.
fluidReactionThermophysicalTransportModel & thermophysicalTransport
Foam::fvModels & fvModels
messageStream Info
Execute application functionObjects to post-process existing results.
Read the control parameters used by setDeltaT.