XiModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Class
25  Foam::XiModel
26 
27 Description
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.
30 
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).
34 
35  See report TR/HGW/10 for details on the Weller two equations model.
36 
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.
40 
41  \f$\Xi_{eq}\f$ is calculated as follows:
42 
43  \f$\Xi_{eq} = 1 + (1 + 2\Xi_{coeff}(0.5 - \dwea{b}))(\Xi^* - 1)\f$
44 
45  where:
46 
47  \f$ \dwea{b} \f$ is the regress variable.
48 
49  \f$ \Xi_{coeff} \f$ is a model constant.
50 
51  \f$ \Xi^* \f$ is the total equilibrium wrinkling combining the effects
52  of the flame inestability and turbulence interaction and is given by
53 
54  \f[
55  \Xi^* = \frac {R}{R - G_\eta - G_{in}}
56  \f]
57 
58  where:
59 
60  \f$ G_\eta \f$ is the generation rate of wrinkling due to turbulence
61  interaction.
62 
63  \f$ G_{in} = \kappa \rho_{u}/\rho_{b} \f$ is the generation
64  rate due to the flame inestability.
65 
66  By adding the removal rates of the two effects:
67 
68  \f[
69  R = G_\eta \frac{\Xi_{\eta_{eq}}}{\Xi_{\eta_{eq}} - 1}
70  + G_{in} \frac{\Xi_{{in}_{eq}}}{\Xi_{{in}_{eq}} - 1}
71  \f]
72 
73  where:
74 
75  \f$ R \f$ is the total removal.
76 
77  \f$ G_\eta \f$ is a model constant.
78 
79  \f$ \Xi_{\eta_{eq}} \f$ is the flame wrinkling due to turbulence.
80 
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).
83 
84 
85 SourceFiles
86  XiModel.C
87 
88 \*---------------------------------------------------------------------------*/
89 
90 #ifndef XiModel_H
91 #define XiModel_H
92 
93 #include "IOdictionary.H"
94 #include "psiuReactionThermo.H"
97 #include "fvcDiv.H"
98 #include "runTimeSelectionTables.H"
99 
100 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
101 
102 namespace Foam
103 {
104 
105 /*---------------------------------------------------------------------------*\
106  Class XiModel Declaration
107 \*---------------------------------------------------------------------------*/
109 class XiModel
110 {
111 
112 protected:
113 
114  // Protected data
123  const surfaceScalarField& phi_;
124 
125  //- Flame wrinking field
127 
128 
129 private:
130 
131  // Private Member Functions
132 
133  //- Disallow copy construct
134  XiModel(const XiModel&);
135 
136  //- Disallow default bitwise assignment
137  void operator=(const XiModel&);
138 
139 
140 public:
141 
142  //- Runtime type information
143  TypeName("XiModel");
144 
145 
146  // Declare run-time constructor selection table
147 
149  (
150  autoPtr,
151  XiModel,
152  dictionary,
153  (
154  const dictionary& XiProperties,
155  const psiuReactionThermo& thermo,
157  const volScalarField& Su,
158  const volScalarField& rho,
159  const volScalarField& b,
160  const surfaceScalarField& phi
161  ),
162  (
163  XiProperties,
164  thermo,
165  turbulence,
166  Su,
167  rho,
168  b,
169  phi
170  )
171  );
172 
173 
174  // Selectors
175 
176  //- Return a reference to the selected Xi model
177  static autoPtr<XiModel> New
178  (
179  const dictionary& XiProperties,
180  const psiuReactionThermo& thermo,
181  const compressible::RASModel& turbulence,
182  const volScalarField& Su,
183  const volScalarField& rho,
184  const volScalarField& b,
185  const surfaceScalarField& phi
186  );
187 
188 
189  // Constructors
190 
191  //- Construct from components
192  XiModel
193  (
194  const dictionary& XiProperties,
195  const psiuReactionThermo& thermo,
196  const compressible::RASModel& turbulence,
197  const volScalarField& Su,
198  const volScalarField& rho,
199  const volScalarField& b,
200  const surfaceScalarField& phi
201  );
202 
203 
204  //- Destructor
205  virtual ~XiModel();
206 
207 
208  // Member Functions
209 
210  //- Return the flame-wrinking Xi
211  virtual const volScalarField& Xi() const
212  {
213  return Xi_;
214  }
215 
216  //- Return the flame diffusivity
217  virtual tmp<volScalarField> Db() const
218  {
219  return turbulence_.muEff();
220  }
221 
222  //- Add Xi to the multivariateSurfaceInterpolationScheme table
223  // if required
224  virtual void addXi
225  (
227  )
228  {}
229 
230  //- Correct the flame-wrinking Xi
231  virtual void correct() = 0;
232 
233  //- Correct the flame-wrinking Xi using the given convection scheme
234  virtual void correct(const fv::convectionScheme<scalar>&)
235  {
236  correct();
237  }
238 
239  //- Update properties from given dictionary
240  virtual bool read(const dictionary& XiProperties) = 0;
241 
242  //- Write fields related to Xi model
243  virtual void writeFields() = 0;
244 };
245 
246 
247 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 
249 } // End namespace Foam
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 #endif
254 
255 // ************************************************************************* //
zeroField Su
Definition: alphaSuSp.H:1
surfaceScalarField & phi
const volScalarField & b_
Definition: XiModel.H:121
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const compressible::RASModel & turbulence_
Definition: XiModel.H:118
Base-class for all Xi models used by the b-Xi combustion model. See Technical Report SH/RE/01R for de...
Definition: XiModel.H:108
virtual ~XiModel()
Destructor.
virtual void writeFields()=0
Write fields related to Xi model.
virtual const volScalarField & Xi() const
Return the flame-wrinking Xi.
Definition: XiModel.H:210
static autoPtr< XiModel > New(const dictionary &XiProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su, const volScalarField &rho, const volScalarField &b, const surfaceScalarField &phi)
Return a reference to the selected Xi model.
volScalarField Xi_
Flame wrinking field.
Definition: XiModel.H:125
Generic GeometricField class.
const surfaceScalarField & phi_
Definition: XiModel.H:122
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
Definition: XiModel.H:216
rhoReactionThermo & thermo
Definition: createFields.H:28
virtual bool read(const dictionary &XiProperties)=0
Update properties from given dictionary.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
const volScalarField & rho_
Definition: XiModel.H:120
declareRunTimeSelectionTable(autoPtr, XiModel, dictionary,(const dictionary &XiProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su, const volScalarField &rho, const volScalarField &b, const surfaceScalarField &phi),(XiProperties, thermo, turbulence, Su, rho, b, phi))
const volScalarField & Su_
Definition: XiModel.H:119
dictionary XiModelCoeffs_
Definition: XiModel.H:115
Foam::psiuReactionThermo.
Calculate the divergence of the given field.
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
Abstract base class for multi-variate surface interpolation schemes.
virtual void correct()=0
Correct the flame-wrinking Xi.
TypeName("XiModel")
Runtime type information.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual void addXi(multivariateSurfaceInterpolationScheme< scalar >::fieldTable &)
Add Xi to the multivariateSurfaceInterpolationScheme table.
Definition: XiModel.H:224
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
const psiuReactionThermo & thermo_
Definition: XiModel.H:117