BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / applications / solvers / combustion / PDRFoam / XiModels / XiModel / XiModel.H
blob8509f05a5c697324f502c043e2c91b6b091c361b
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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/>.
24 Class
25     Foam::XiModel
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.
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$
45     where:
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
54         \f[
55             \Xi^* = \frac {R}{R - G_\eta - G_{in}}
56         \f]
58     where:
60         \f$ G_\eta \f$ is the generation rate of wrinkling due to turbulence
61         interaction.
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:
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]
73     where:
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).
85 SourceFiles
86     XiModel.C
88 \*---------------------------------------------------------------------------*/
90 #ifndef XiModel_H
91 #define XiModel_H
93 #include "IOdictionary.H"
94 #include "hhuCombustionThermo.H"
95 #include "RASModel.H"
96 #include "multivariateSurfaceInterpolationScheme.H"
97 #include "runTimeSelectionTables.H"
99 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
101 namespace Foam
104 /*---------------------------------------------------------------------------*\
105                           Class XiModel Declaration
106 \*---------------------------------------------------------------------------*/
108 class XiModel
111 protected:
113     // Protected data
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
125         volScalarField Xi_;
128 private:
130     // Private Member Functions
132         //- Disallow copy construct
133         XiModel(const XiModel&);
135         //- Disallow default bitwise assignment
136         void operator=(const XiModel&);
139 public:
141     //- Runtime type information
142     TypeName("XiModel");
145     // Declare run-time constructor selection table
147         declareRunTimeSelectionTable
148         (
149             autoPtr,
150             XiModel,
151             dictionary,
152             (
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
160             ),
161             (
162                 XiProperties,
163                 thermo,
164                 turbulence,
165                 Su,
166                 rho,
167                 b,
168                 phi
169             )
170         );
173     // Selectors
175         //- Return a reference to the selected Xi model
176         static autoPtr<XiModel> New
177         (
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
185         );
188     // Constructors
190         //- Construct from components
191         XiModel
192         (
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
200         );
203     //- Destructor
204     virtual ~XiModel();
207     // Member Functions
209         //- Return the flame-wrinking Xi
210         virtual const volScalarField& Xi() const
211         {
212             return Xi_;
213         }
215         //- Return the flame diffusivity
216         virtual tmp<volScalarField> Db() const
217         {
218             return turbulence_.muEff();
219         }
221         //- Add Xi to the multivariateSurfaceInterpolationScheme table
222         //  if required
223         virtual void addXi
224         (
225             multivariateSurfaceInterpolationScheme<scalar>::fieldTable&
226         )
227         {}
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>&)
234         {
235             correct();
236         }
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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 #endif
254 // ************************************************************************* //