1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
29 Compressible premixed/partially-premixed combustion solver with turbulence
32 Combusting RANS code using the b-Xi two-equation model.
33 Xi may be obtained by either the solution of the Xi transport
34 equation or from an algebraic exression. Both approaches are
35 based on Gulder's flame speed correlation which has been shown
36 to be appropriate by comparison with the results from the
39 Strain effects are encorporated directly into the Xi equation
40 but not in the algebraic approximation. Further work need to be
41 done on this issue, particularly regarding the enhanced removal rate
42 caused by flame compression. Analysis using results of the spectral
43 model will be required.
45 For cases involving very lean Propane flames or other flames which are
46 very strain-sensitive, a transport equation for the laminar flame
47 speed is present. This equation is derived using heuristic arguments
48 involving the strain time scale and the strain-rate at extinction.
49 the transport velocity is the same as that for the Xi equation.
51 For large flames e.g. explosions additional modelling for the flame
52 wrinkling due to surface instabilities may be applied.
54 PDR (porosity/distributed resistance) modelling is included to handle
55 regions containing blockages which cannot be resolved by the mesh.
57 \*---------------------------------------------------------------------------*/
59 #include <finiteVolume/fvCFD.H>
60 #include <dynamicFvMesh/dynamicFvMesh.H>
61 #include <reactionThermophysicalModels/hhuCombustionThermo.H>
62 #include <compressibleRASModels/RASModel.H>
63 #include <laminarFlameSpeedModels/laminarFlameSpeed.H>
65 #include "PDRDragModel.H"
66 #include <engine/ignition.H>
67 #include <OpenFOAM/Switch.H>
68 #include <finiteVolume/bound.H>
69 #include <dynamicFvMesh/dynamicRefineFvMesh.H>
71 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 int main(int argc, char *argv[])
75 #include <OpenFOAM/setRootCase.H>
77 #include <OpenFOAM/createTime.H>
78 #include <dynamicFvMesh/createDynamicFvMesh.H>
79 #include "readCombustionProperties.H"
80 #include <finiteVolume/readGravitationalAcceleration.H>
81 #include "createFields.H"
82 #include <finiteVolume/readPISOControls.H>
83 #include <finiteVolume/initContinuityErrs.H>
84 #include <finiteVolume/readTimeControls.H>
85 #include <finiteVolume/setInitialDeltaT.H>
89 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91 Info<< "\nStarting time loop\n" << endl;
95 #include <finiteVolume/readTimeControls.H>
96 #include <finiteVolume/readPISOControls.H>
97 #include <finiteVolume/CourantNo.H>
99 #include "setDeltaT.H"
103 Info<< "\n\nTime = " << runTime.timeName() << endl;
105 // Indicators for refinement. Note: before runTime++
106 // only for postprocessing reasons.
107 tmp<volScalarField> tmagGradP = mag(fvc::grad(p));
108 volScalarField normalisedGradP
111 tmagGradP()/max(tmagGradP())
113 normalisedGradP.writeOpt() = IOobject::AUTO_WRITE;
116 bool meshChanged = false;
118 // Make the fluxes absolute
119 fvc::makeAbsolute(phi, rho, U);
121 // Test : disable refinement for some cells
122 PackedBoolList& protectedCell =
123 refCast<dynamicRefineFvMesh>(mesh).protectedCell();
125 if (protectedCell.empty())
127 protectedCell.setSize(mesh.nCells());
133 if (betav[cellI] < 0.99)
135 protectedCell[cellI] = 1;
139 //volScalarField pIndicator("pIndicator",
140 // p*(fvc::laplacian(p))
142 // magSqr(fvc::grad(p))
143 // + dimensionedScalar
146 // sqr(p.dimensions()/dimLength),
150 //pIndicator.writeOpt() = IOobject::AUTO_WRITE;
152 // Flux estimate for introduced faces.
153 volVectorField rhoU("rhoU", rho*U);
155 // Do any mesh changes
156 meshChanged = mesh.update();
158 // if (mesh.moving() || meshChanged)
160 // #include "correctPhi.H"
163 // Make the fluxes relative to the mesh motion
164 fvc::makeRelative(phi, rho, U);
172 for (int corr=1; corr<=nCorr; corr++)
187 turbulence->correct();
191 Info<< "\nExecutionTime = "
192 << runTime.elapsedCpuTime()
202 // ************************ vim: set sw=4 sts=4 et: ************************ //