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-2024 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 
30  References:
31  \verbatim
32  Weller, H. G. (1993).
33  The development of a new flame area combustion model
34  using conditional averaging.
35  Thermo-fluids section report TF 9307.
36 
37  Weller, H. G., Tabor, G., Gosman, A. D., & Fureby, C. (1998, January).
38  Application of a flame-wrinkling LES combustion model
39  to a turbulent mixing layer.
40  In Symposium (International) on combustion
41  (Vol. 27, No. 1, pp. 899-907). Elsevier.
42  \endverbatim
43 
44  Xi is given through an algebraic expression (Foam::XiModels::algebraic), by
45  solving a transport equation (Foam::XiModels::transport.H) or a fixed value
46  (Foam::XiModels::fixed).
47 
48  In the algebraic and transport methods \f$\Xi_{eq}\f$ is calculated in
49  similar way. In the algebraic approach, \f$\Xi_{eq}\f$ is the value used in
50  the \f$ b \f$ transport equation.
51 
52  \f$\Xi_{eq}\f$ is calculated as follows:
53 
54  \f$\Xi_{eq} = 1 + (1 + 2\Xi_{coeff}(0.5 - \dwea{b}))(\Xi^* - 1)\f$
55 
56  where:
57 
58  \f$ \dwea{b} \f$ is the regress variable.
59 
60  \f$ \Xi_{coeff} \f$ is a model constant.
61 
62  \f$ \Xi^* \f$ is the total equilibrium wrinkling combining the effects
63  of the flame instability and turbulence interaction and is given by
64 
65  \f[
66  \Xi^* = \frac {R}{R - G_\eta - G_{in}}
67  \f]
68 
69  where:
70 
71  \f$ G_\eta \f$ is the generation rate of wrinkling due to turbulence
72  interaction.
73 
74  \f$ G_{in} = \kappa \rho_{u}/\rho_{b} \f$ is the generation
75  rate due to the flame instability.
76 
77  By adding the removal rates of the two effects:
78 
79  \f[
80  R = G_\eta \frac{\Xi_{\eta_{eq}}}{\Xi_{\eta_{eq}} - 1}
81  + G_{in} \frac{\Xi_{{in}_{eq}}}{\Xi_{{in}_{eq}} - 1}
82  \f]
83 
84  where:
85 
86  \f$ R \f$ is the total removal.
87 
88  \f$ G_\eta \f$ is a model constant.
89 
90  \f$ \Xi_{\eta_{eq}} \f$ is the flame wrinkling due to turbulence.
91 
92  \f$ \Xi_{{in}_{eq}} \f$ is the equilibrium level of the flame wrinkling
93  generated by instability. It is a constant (default 2.5).
94 
95 
96 SourceFiles
97  XiModel.C
98 
99 \*---------------------------------------------------------------------------*/
100 
101 #ifndef XiModel_H
102 #define XiModel_H
103 
106 #include "runTimeSelectionTables.H"
107 
108 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109 
110 namespace Foam
111 {
112 
113 /*---------------------------------------------------------------------------*\
114  Class XiModel Declaration
115 \*---------------------------------------------------------------------------*/
116 
117 class XiModel
118 {
119 
120 protected:
121 
122  // Protected data
123 
124  //- Thermo
126 
127  //- Thermo-physical transport
129 
130  //- Turbulence
132 
133  const volScalarField& Su_;
134 
135  const volScalarField& rho_;
136 
137  const volScalarField& b_;
138 
139  //- Flame wrinkling field
141 
142 
143  // Protected member functions
144 
145  //- Update coefficients from given dictionary
146  virtual bool readCoeffs(const dictionary& dict) = 0;
147 
148 
149 public:
150 
151  //- Runtime type information
152  TypeName("XiModel");
153 
154 
155  // Declare run-time constructor selection table
156 
158  (
159  autoPtr,
160  XiModel,
161  dictionary,
162  (
163  const dictionary& dict,
166  const volScalarField& Su
167  ),
168  (
169  dict,
170  thermo,
171  turbulence,
172  Su
173  )
174  );
175 
176 
177  // Constructors
178 
179  //- Construct from components
180  XiModel
181  (
184  const volScalarField& Su
185  );
186 
187  //- Disallow default bitwise copy construction
188  XiModel(const XiModel&) = delete;
189 
190 
191  // Selectors
192 
193  //- Return a reference to the selected Xi model
194  static autoPtr<XiModel> New
195  (
196  const dictionary& combustionProperties,
199  const volScalarField& Su
200  );
201 
202 
203  //- Destructor
204  virtual ~XiModel();
205 
206 
207  // Member Functions
208 
209  //- Return the flame-wrinkling Xi
210  virtual const volScalarField& Xi() const
211  {
212  return Xi_;
213  }
214 
215  //- Return the flame diffusivity
216  virtual tmp<volScalarField> Db() const
217  {
218  return volScalarField::New
219  (
220  "Db",
222  );
223  }
224 
225  //- Add Xi to the multivariateSurfaceInterpolationScheme table
226  // if required
227  virtual void addXi
228  (
230  )
231  {}
232 
233  //- Reset Xi to the unburnt state (1)
234  virtual void reset();
235 
236  //- Correct the flame-wrinkling Xi
237  virtual void correct() = 0;
238 
239  //- Update properties from the given combustionProperties dictionary
240  bool read(const dictionary& combustionProperties);
241 
242 
243  // Member Operators
244 
245  //- Disallow default bitwise assignment
246  void operator=(const XiModel&) = delete;
247 };
248 
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 } // End namespace Foam
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 #endif
257 
258 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Templated abstract base class for thermophysical transport models.
Base-class for all Xi models used by the b-Xi combustion model.
Definition: XiModel.H:117
virtual ~XiModel()
Destructor.
Definition: XiModel.C:78
const volScalarField & rho_
Definition: XiModel.H:134
const volScalarField & Su_
Definition: XiModel.H:132
void operator=(const XiModel &)=delete
Disallow default bitwise assignment.
bool read(const dictionary &combustionProperties)
Update properties from the given combustionProperties dictionary.
Definition: XiModel.C:84
const psiuMulticomponentThermo & thermo_
Thermo.
Definition: XiModel.H:124
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
Definition: XiModel.H:215
virtual void correct()=0
Correct the flame-wrinkling Xi.
const volScalarField & b_
Definition: XiModel.H:136
const compressibleMomentumTransportModel & turbulence_
Turbulence.
Definition: XiModel.H:130
virtual const volScalarField & Xi() const
Return the flame-wrinkling Xi.
Definition: XiModel.H:209
volScalarField Xi_
Flame wrinkling field.
Definition: XiModel.H:139
declareRunTimeSelectionTable(autoPtr, XiModel, dictionary,(const dictionary &dict, const psiuMulticomponentThermo &thermo, const fluidThermoThermophysicalTransportModel &turbulence, const volScalarField &Su),(dict, thermo, turbulence, Su))
TypeName("XiModel")
Runtime type information.
const fluidThermoThermophysicalTransportModel & thermoTransport_
Thermo-physical transport.
Definition: XiModel.H:127
virtual void addXi(multivariateSurfaceInterpolationScheme< scalar >::fieldTable &)
Add Xi to the multivariateSurfaceInterpolationScheme table.
Definition: XiModel.H:227
virtual void reset()
Reset Xi to the unburnt state (1)
Definition: XiModel.C:93
virtual bool readCoeffs(const dictionary &dict)=0
Update coefficients from given dictionary.
Definition: XiModel.C:39
static autoPtr< XiModel > New(const dictionary &combustionProperties, const psiuMulticomponentThermo &thermo, const fluidThermoThermophysicalTransportModel &turbulence, const volScalarField &Su)
Return a reference to the selected Xi model.
Definition: XiModelNew.C:31
XiModel(const psiuMulticomponentThermo &thermo, const fluidThermoThermophysicalTransportModel &turbulence, const volScalarField &Su)
Construct from components.
Definition: XiModel.C:48
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Base class for single-phase compressible turbulence models.
const rhoField & rho() const
Return the density field.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
Abstract base class for multi-variate surface interpolation schemes.
Base-class for combustion fluid thermodynamic properties based on compressibility.
A class for managing temporary objects.
Definition: tmp.H:55
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
Namespace for OpenFOAM.
Macros to ease declaration of run-time selection tables.
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
fluidMulticomponentThermo & thermo
Definition: createFields.H:31