Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / applications / solvers / combustion / PDRFoam / PDRFoam.C
blob3c4da27187bf089e0a534ab4738b3d0d8dc05941
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2011 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 Application
25     PDRFoam
27 Description
28     Solver for compressible premixed/partially-premixed combustion with
29     turbulence modelling.
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.
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.
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.
50     For large flames e.g. explosions additional modelling for the flame
51     wrinkling due to surface instabilities may be applied.
53     PDR (porosity/distributed resistance) modelling is included to handle
54     regions containing blockages which cannot be resolved by the mesh.
56     The fields used by this solver are:
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.
68 \*---------------------------------------------------------------------------*/
70 #include "fvCFD.H"
71 #include "hhuCombustionThermo.H"
72 #include "RASModel.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"
81 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83 int main(int argc, char *argv[])
85     #include "setRootCase.H"
87     #include "createTime.H"
88     #include "createMesh.H"
89     #include "readCombustionProperties.H"
90     #include "readGravitationalAcceleration.H"
91     #include "createFields.H"
92     #include "initContinuityErrs.H"
93     #include "readTimeControls.H"
94     #include "compressibleCourantNo.H"
95     #include "setInitialDeltaT.H"
97     pimpleControl pimple(mesh);
99     scalar StCoNum = 0.0;
101     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
103     Info<< "\nStarting time loop\n" << endl;
105     while (runTime.run())
106     {
107         #include "readTimeControls.H"
108         #include "compressibleCourantNo.H"
109         #include "setDeltaT.H"
111         runTime++;
112         Info<< "\n\nTime = " << runTime.timeName() << endl;
114         #include "rhoEqn.H"
116         // --- Pressure-velocity PIMPLE corrector loop
117         for (pimple.start(); pimple.loop(); pimple++)
118         {
119             #include "UEqn.H"
121             // --- PISO loop
122             for (int corr=1; corr<=pimple.nCorr(); corr++)
123             {
124                 #include "bEqn.H"
125                 #include "ftEqn.H"
126                 #include "huEqn.H"
127                 #include "hEqn.H"
129                 if (!ign.ignited())
130                 {
131                     hu == h;
132                 }
134                 #include "pEqn.H"
135             }
137             if (pimple.turbCorr())
138             {
139                 turbulence->correct();
140             }
141         }
143         runTime.write();
145         Info<< "\nExecutionTime = "
146              << runTime.elapsedCpuTime()
147              << " s\n" << endl;
148     }
150     Info<< "\n end\n";
152     return 0;
156 // ************************************************************************* //