1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
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
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
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/>.
28 Solver for compressible premixed/partially-premixed combustion with
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
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 \*---------------------------------------------------------------------------*/
59 #include "dynamicFvMesh.H"
60 #include "hhuCombustionThermo.H"
62 #include "laminarFlameSpeed.H"
64 #include "PDRDragModel.H"
68 #include "dynamicRefineFvMesh.H"
69 #include "pimpleControl.H"
71 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 int main(int argc, char *argv[])
75 #include "setRootCase.H"
77 #include "createTime.H"
78 #include "createDynamicFvMesh.H"
79 #include "readCombustionProperties.H"
80 #include "readGravitationalAcceleration.H"
81 #include "createFields.H"
82 #include "initContinuityErrs.H"
83 #include "readTimeControls.H"
84 #include "compressibleCourantNo.H"
85 #include "setInitialDeltaT.H"
87 pimpleControl pimple(mesh);
91 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
93 Info<< "\nStarting time loop\n" << endl;
95 bool hasChanged = false;
99 #include "readTimeControls.H"
100 #include "compressibleCourantNo.H"
101 #include "setDeltaT.H"
103 // Indicators for refinement. Note: before runTime++
104 // only for postprocessing reasons.
105 tmp<volScalarField> tmagGradP = mag(fvc::grad(p));
106 volScalarField normalisedGradP
109 tmagGradP()/max(tmagGradP())
111 normalisedGradP.writeOpt() = IOobject::AUTO_WRITE;
116 Info<< "\n\nTime = " << runTime.timeName() << endl;
119 // Make the fluxes absolute
120 fvc::makeAbsolute(phi, rho, U);
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 // Flux estimate for introduced faces.
140 volVectorField rhoU("rhoU", rho*U);
142 // Do any mesh changes
143 bool meshChanged = mesh.update();
150 if (runTime.write() && hasChanged)
156 flameWrinkling->writeFields();
160 // Make the fluxes relative to the mesh motion
161 fvc::makeRelative(phi, rho, U);
167 // --- Pressure-velocity PIMPLE corrector loop
168 for (pimple.start(); pimple.loop(); pimple++)
174 for (int corr=1; corr<=pimple.nCorr(); corr++)
189 if (pimple.turbCorr())
191 turbulence->correct();
197 Info<< "\nExecutionTime = "
198 << runTime.elapsedCpuTime()
208 // ************************************************************************* //