00001 /*---------------------------------------------------------------------------*\ 00002 ========= | 00003 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox 00004 \\ / O peration | 00005 \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd. 00006 \\/ M anipulation | 00007 ------------------------------------------------------------------------------- 00008 License 00009 This file is part of OpenFOAM. 00010 00011 OpenFOAM is free software: you can redistribute it and/or modify it 00012 under the terms of the GNU General Public License as published by 00013 the Free Software Foundation, either version 3 of the License, or 00014 (at your option) any later version. 00015 00016 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT 00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or 00018 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License 00019 for more details. 00020 00021 You should have received a copy of the GNU General Public License 00022 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>. 00023 00024 Class 00025 Foam::XiModel 00026 00027 Description 00028 Base-class for all Xi models used by the b-Xi combustion model. 00029 See Technical Report SH/RE/01R for details on the PDR modelling. 00030 00031 Xi is given through an algebraic expression (\link algebraic.H \endlink), 00032 by solving a transport equation (\link transport.H \endlink) or a 00033 fixed value (\link fixed.H \endlink). 00034 00035 See report TR/HGW/10 for details on the Weller two equations model. 00036 00037 In the algebraic and transport methods \f$\Xi_{eq}\f$ is calculated in 00038 similar way. In the algebraic approach, \f$\Xi_{eq}\f$ is the value used in 00039 the \f$ b \f$ transport equation. 00040 00041 \f$\Xi_{eq}\f$ is calculated as follows: 00042 00043 \f$\Xi_{eq} = 1 + (1 + 2\Xi_{coeff}(0.5 - \dwea{b}))(\Xi^* - 1)\f$ 00044 00045 where: 00046 00047 \f$ \dwea{b} \f$ is the regress variable. 00048 00049 \f$ \Xi_{coeff} \f$ is a model constant. 00050 00051 \f$ \Xi^* \f$ is the total equilibrium wrinkling combining the effects 00052 of the flame inestability and turbulence interaction and is given by 00053 00054 \f[ 00055 \Xi^* = \frac {R}{R - G_\eta - G_{in}} 00056 \f] 00057 00058 where: 00059 00060 \f$ G_\eta \f$ is the generation rate of wrinkling due to turbulence 00061 interaction. 00062 00063 \f$ G_{in} = \kappa \rho_{u}/\rho_{b} \f$ is the generation 00064 rate due to the flame inestability. 00065 00066 By adding the removal rates of the two effects: 00067 00068 \f[ 00069 R = G_\eta \frac{\Xi_{\eta_{eq}}}{\Xi_{\eta_{eq}} - 1} 00070 + G_{in} \frac{\Xi_{{in}_{eq}}}{\Xi_{{in}_{eq}} - 1} 00071 \f] 00072 00073 where: 00074 00075 \f$ R \f$ is the total removal. 00076 00077 \f$ G_\eta \f$ is a model constant. 00078 00079 \f$ \Xi_{\eta_{eq}} \f$ is the flame wrinkling due to turbulence. 00080 00081 \f$ \Xi_{{in}_{eq}} \f$ is the equilibrium level of the flame wrinkling 00082 generated by inestability. It is a constant (default 2.5). 00083 00084 00085 SourceFiles 00086 XiModel.C 00087 00088 \*---------------------------------------------------------------------------*/ 00089 00090 #ifndef XiModel_H 00091 #define XiModel_H 00092 00093 #include <OpenFOAM/IOdictionary.H> 00094 #include <reactionThermophysicalModels/hhuCombustionThermo.H> 00095 #include <compressibleRASModels/RASModel.H> 00096 #include <finiteVolume/multivariateSurfaceInterpolationScheme.H> 00097 #include <OpenFOAM/runTimeSelectionTables.H> 00098 00099 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00100 00101 namespace Foam 00102 { 00103 00104 /*---------------------------------------------------------------------------*\ 00105 Class XiModel Declaration 00106 \*---------------------------------------------------------------------------*/ 00107 00108 class XiModel 00109 { 00110 00111 protected: 00112 00113 // Protected data 00114 00115 dictionary XiModelCoeffs_; 00116 00117 const hhuCombustionThermo& thermo_; 00118 const compressible::RASModel& turbulence_; 00119 const volScalarField& Su_; 00120 const volScalarField& rho_; 00121 const volScalarField& b_; 00122 const surfaceScalarField& phi_; 00123 00124 //- Flame wrinking field 00125 volScalarField Xi_; 00126 00127 00128 private: 00129 00130 // Private Member Functions 00131 00132 //- Disallow copy construct 00133 XiModel(const XiModel&); 00134 00135 //- Disallow default bitwise assignment 00136 void operator=(const XiModel&); 00137 00138 00139 public: 00140 00141 //- Runtime type information 00142 TypeName("XiModel"); 00143 00144 00145 // Declare run-time constructor selection table 00146 00147 declareRunTimeSelectionTable 00148 ( 00149 autoPtr, 00150 XiModel, 00151 dictionary, 00152 ( 00153 const dictionary& XiProperties, 00154 const hhuCombustionThermo& thermo, 00155 const compressible::RASModel& turbulence, 00156 const volScalarField& Su, 00157 const volScalarField& rho, 00158 const volScalarField& b, 00159 const surfaceScalarField& phi 00160 ), 00161 ( 00162 XiProperties, 00163 thermo, 00164 turbulence, 00165 Su, 00166 rho, 00167 b, 00168 phi 00169 ) 00170 ); 00171 00172 00173 // Selectors 00174 00175 //- Return a reference to the selected Xi model 00176 static autoPtr<XiModel> New 00177 ( 00178 const dictionary& XiProperties, 00179 const hhuCombustionThermo& thermo, 00180 const compressible::RASModel& turbulence, 00181 const volScalarField& Su, 00182 const volScalarField& rho, 00183 const volScalarField& b, 00184 const surfaceScalarField& phi 00185 ); 00186 00187 00188 // Constructors 00189 00190 //- Construct from components 00191 XiModel 00192 ( 00193 const dictionary& XiProperties, 00194 const hhuCombustionThermo& thermo, 00195 const compressible::RASModel& turbulence, 00196 const volScalarField& Su, 00197 const volScalarField& rho, 00198 const volScalarField& b, 00199 const surfaceScalarField& phi 00200 ); 00201 00202 00203 // Destructor 00204 00205 virtual ~XiModel(); 00206 00207 00208 // Member Functions 00209 00210 //- Return the flame-wrinking Xi 00211 virtual const volScalarField& Xi() const 00212 { 00213 return Xi_; 00214 } 00215 00216 //- Return the flame diffusivity 00217 virtual tmp<volScalarField> Db() const 00218 { 00219 return turbulence_.muEff(); 00220 } 00221 00222 //- Add Xi to the multivariateSurfaceInterpolationScheme table 00223 // if required 00224 virtual void addXi 00225 ( 00226 multivariateSurfaceInterpolationScheme<scalar>::fieldTable& 00227 ) 00228 {} 00229 00230 //- Correct the flame-wrinking Xi 00231 virtual void correct() = 0; 00232 00233 //- Correct the flame-wrinking Xi using the given convection scheme 00234 virtual void correct(const fv::convectionScheme<scalar>&) 00235 { 00236 correct(); 00237 } 00238 00239 //- Update properties from given dictionary 00240 virtual bool read(const dictionary& XiProperties) = 0; 00241 }; 00242 00243 00244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00245 00246 } // End namespace Foam 00247 00248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00249 00250 #endif 00251 00252 // ************************ vim: set sw=4 sts=4 et: ************************ //