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-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::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 private:
70 
71  // Private Member Functions
72 
73  //- Disallow copy construct
74  XiGModel(const XiGModel&);
75 
76  //- Disallow default bitwise assignment
77  void operator=(const XiGModel&);
78 
79 
80 public:
81 
82  //- Runtime type information
83  TypeName("XiGModel");
84 
85 
86  // Declare run-time constructor selection table
87 
89  (
90  autoPtr,
91  XiGModel,
92  dictionary,
93  (
94  const dictionary& XiGProperties,
97  const volScalarField& Su
98  ),
99  (
100  XiGProperties,
101  thermo,
102  turbulence,
103  Su
104  )
105  );
106 
107 
108  // Selectors
109 
110  //- Return a reference to the selected XiG model
111  static autoPtr<XiGModel> New
112  (
113  const dictionary& XiGProperties,
114  const psiuReactionThermo& thermo,
115  const compressible::RASModel& turbulence,
116  const volScalarField& Su
117  );
118 
119 
120  // Constructors
121 
122  //- Construct from components
123  XiGModel
124  (
125  const dictionary& XiGProperties,
126  const psiuReactionThermo& thermo,
127  const compressible::RASModel& turbulence,
128  const volScalarField& Su
129  );
130 
131 
132  //- Destructor
133  virtual ~XiGModel();
134 
135 
136  // Member Functions
137 
138  //- Return the flame-wrinking genration rate
139  virtual tmp<volScalarField> G() const = 0;
140 
141  //- Return the flame diffusivity
142  virtual tmp<volScalarField> Db() const
143  {
144  return turbulence_.muEff();
145  }
146 
147  //- Update properties from given dictionary
148  virtual bool read(const dictionary& XiGProperties) = 0;
149 };
150 
151 
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 
154 } // End namespace Foam
155 
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 
158 #endif
159 
160 // ************************************************************************* //
virtual ~XiGModel()
Destructor.
zeroField Su
Definition: alphaSuSp.H:1
const psiuReactionThermo & thermo_
Definition: XiGModel.H:63
virtual tmp< volScalarField > G() const =0
Return the flame-wrinking genration rate.
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: XiGModel.H:64
rhoReactionThermo & thermo
Definition: createFields.H:28
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
Foam::psiuReactionThermo.
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("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
virtual bool read(const dictionary &XiGProperties)=0
Update properties from given dictionary.
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:141
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
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.
Namespace for OpenFOAM.