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 Base-class for all Xi models used by the b-Xi combustion model.
29 See Technical Report SH/RE/01R for details on the PDR modelling.
31 Xi is given through an algebraic expression (\link algebraic.H \endlink),
32 by solving a transport equation (\link transport.H \endlink) or a
33 fixed value (\link fixed.H \endlink).
35 See report TR/HGW/10 for details on the Weller two equations model.
37 In the algebraic and transport methods \f$\Xi_{eq}\f$ is calculated in
38 similar way. In the algebraic approach, \f$\Xi_{eq}\f$ is the value used in
39 the \f$ b \f$ transport equation.
41 \f$\Xi_{eq}\f$ is calculated as follows:
43 \f$\Xi_{eq} = 1 + (1 + 2\Xi_{coeff}(0.5 - \dwea{b}))(\Xi^* - 1)\f$
47 \f$ \dwea{b} \f$ is the regress variable.
49 \f$ \Xi_{coeff} \f$ is a model constant.
51 \f$ \Xi^* \f$ is the total equilibrium wrinkling combining the effects
52 of the flame inestability and turbulence interaction and is given by
55 \Xi^* = \frac {R}{R - G_\eta - G_{in}}
60 \f$ G_\eta \f$ is the generation rate of wrinkling due to turbulence
63 \f$ G_{in} = \kappa \rho_{u}/\rho_{b} \f$ is the generation
64 rate due to the flame inestability.
66 By adding the removal rates of the two effects:
69 R = G_\eta \frac{\Xi_{\eta_{eq}}}{\Xi_{\eta_{eq}} - 1}
70 + G_{in} \frac{\Xi_{{in}_{eq}}}{\Xi_{{in}_{eq}} - 1}
75 \f$ R \f$ is the total removal.
77 \f$ G_\eta \f$ is a model constant.
79 \f$ \Xi_{\eta_{eq}} \f$ is the flame wrinkling due to turbulence.
81 \f$ \Xi_{{in}_{eq}} \f$ is the equilibrium level of the flame wrinkling
82 generated by inestability. It is a constant (default 2.5).
88 \*---------------------------------------------------------------------------*/
93 #include "IOdictionary.H"
94 #include "hhuCombustionThermo.H"
96 #include "multivariateSurfaceInterpolationScheme.H"
97 #include "runTimeSelectionTables.H"
99 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104 /*---------------------------------------------------------------------------*\
105 Class XiModel Declaration
106 \*---------------------------------------------------------------------------*/
115 dictionary XiModelCoeffs_;
117 const hhuCombustionThermo& thermo_;
118 const compressible::RASModel& turbulence_;
119 const volScalarField& Su_;
120 const volScalarField& rho_;
121 const volScalarField& b_;
122 const surfaceScalarField& phi_;
124 //- Flame wrinking field
130 // Private Member Functions
132 //- Disallow copy construct
133 XiModel(const XiModel&);
135 //- Disallow default bitwise assignment
136 void operator=(const XiModel&);
141 //- Runtime type information
145 // Declare run-time constructor selection table
147 declareRunTimeSelectionTable
153 const dictionary& XiProperties,
154 const hhuCombustionThermo& thermo,
155 const compressible::RASModel& turbulence,
156 const volScalarField& Su,
157 const volScalarField& rho,
158 const volScalarField& b,
159 const surfaceScalarField& phi
175 //- Return a reference to the selected Xi model
176 static autoPtr<XiModel> New
178 const dictionary& XiProperties,
179 const hhuCombustionThermo& thermo,
180 const compressible::RASModel& turbulence,
181 const volScalarField& Su,
182 const volScalarField& rho,
183 const volScalarField& b,
184 const surfaceScalarField& phi
190 //- Construct from components
193 const dictionary& XiProperties,
194 const hhuCombustionThermo& thermo,
195 const compressible::RASModel& turbulence,
196 const volScalarField& Su,
197 const volScalarField& rho,
198 const volScalarField& b,
199 const surfaceScalarField& phi
209 //- Return the flame-wrinking Xi
210 virtual const volScalarField& Xi() const
215 //- Return the flame diffusivity
216 virtual tmp<volScalarField> Db() const
218 return turbulence_.muEff();
221 //- Add Xi to the multivariateSurfaceInterpolationScheme table
225 multivariateSurfaceInterpolationScheme<scalar>::fieldTable&
229 //- Correct the flame-wrinking Xi
230 virtual void correct() = 0;
232 //- Correct the flame-wrinking Xi using the given convection scheme
233 virtual void correct(const fv::convectionScheme<scalar>&)
238 //- Update properties from given dictionary
239 virtual bool read(const dictionary& XiProperties) = 0;
241 //- Write fields related to Xi model
242 virtual void writeFields() = 0;
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 } // End namespace Foam
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 // ************************************************************************* //