All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
XiGModel.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-2021 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::XiGModel
26 
27 Description
28  Base-class for all Xi generation models used by the b-Xi combustion model.
29  See Technical Report SH/RE/01R for details on the PDR modelling. For details
30  on the use of XiGModel see \link XiModel.H \endlink. The model available is
31  \link instabilityG.H \endlink
32 
33 SourceFiles
34  XiGModel.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef XiGModel_H
39 #define XiGModel_H
40 
41 #include "IOdictionary.H"
42 #include "psiuReactionThermo.H"
44 #include "runTimeSelectionTables.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class XiGModel Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 class XiGModel
56 {
57 
58 protected:
59 
60  // Protected data
61 
63 
66  const volScalarField& Su_;
67 
68 
69 public:
70 
71  //- Runtime type information
72  TypeName("XiGModel");
73 
74 
75  // Declare run-time constructor selection table
76 
78  (
79  autoPtr,
80  XiGModel,
81  dictionary,
82  (
83  const dictionary& XiGProperties,
86  const volScalarField& Su
87  ),
88  (
89  XiGProperties,
90  thermo,
91  turbulence,
92  Su
93  )
94  );
95 
96 
97  // Constructors
98 
99  //- Construct from components
100  XiGModel
101  (
102  const dictionary& XiGProperties,
103  const psiuReactionThermo& thermo,
104  const compressible::RASModel& turbulence,
105  const volScalarField& Su
106  );
107 
108  //- Disallow default bitwise copy construction
109  XiGModel(const XiGModel&);
110 
111 
112  // Selectors
113 
114  //- Return a reference to the selected XiG model
115  static autoPtr<XiGModel> New
116  (
117  const dictionary& XiGProperties,
118  const psiuReactionThermo& thermo,
119  const compressible::RASModel& turbulence,
120  const volScalarField& Su
121  );
122 
123 
124  //- Destructor
125  virtual ~XiGModel();
126 
127 
128  // Member Functions
129 
130  //- Return the flame-wrinkling generation rate
131  virtual tmp<volScalarField> G() const = 0;
132 
133  //- Return the flame diffusivity
134  virtual tmp<volScalarField> Db() const
135  {
136  return turbulence_.muEff();
137  }
138 
139  //- Update properties from given dictionary
140  virtual bool read(const dictionary& XiGProperties) = 0;
141 
142 
143  // Member Operators
144 
145  //- Disallow default bitwise assignment
146  void operator=(const XiGModel&) = delete;
147 };
148 
149 
150 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
151 
152 } // End namespace Foam
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
156 #endif
157 
158 // ************************************************************************* //
virtual ~XiGModel()
Destructor.
fluidReactionThermo & thermo
Definition: createFields.H:28
const psiuReactionThermo & thermo_
Definition: XiGModel.H:63
virtual tmp< volScalarField > G() const =0
Return the flame-wrinkling generation rate.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const compressible::RASModel & turbulence_
Definition: XiGModel.H:64
dictionary XiGModelCoeffs_
Definition: XiGModel.H:61
declareRunTimeSelectionTable(autoPtr, XiGModel, dictionary,(const dictionary &XiGProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su),(XiGProperties, thermo, turbulence, Su))
Base-class for all Xi generation models used by the b-Xi combustion model. See Technical Report SH/RE...
Definition: XiGModel.H:54
Base-class for combustion fluid thermodynamic properties based on compressibility.
const volScalarField & Su_
Definition: XiGModel.H:65
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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
virtual bool read(const dictionary &XiGProperties)=0
Update properties from given dictionary.
XiGModel(const dictionary &XiGProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su)
Construct from components.
TypeName("XiGModel")
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 tmp< volScalarField > Db() const
Return the flame diffusivity.
Definition: XiGModel.H:133
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
static autoPtr< XiGModel > New(const dictionary &XiGProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su)
Return a reference to the selected XiG model.
void operator=(const XiGModel &)=delete
Disallow default bitwise assignment.
Namespace for OpenFOAM.