XiEqModel.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-2020 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::XiEqModel
26 
27 Description
28  Base-class for all XiEq models used by the b-XiEq combustion model.
29  The available models are :
30  \link basicXiSubXiEq.H \endlink
31  \link Gulder.H \endlink
32  \link instabilityXiEq.H \endlink
33  \link SCOPEBlendXiEq.H \endlink
34  \link SCOPEXiEq.H \endlink
35 
36 SourceFiles
37  XiEqModel.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef XiEqModel_H
42 #define XiEqModel_H
43 
44 #include "IOdictionary.H"
45 #include "psiuReactionThermo.H"
47 #include "runTimeSelectionTables.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 /*---------------------------------------------------------------------------*\
55  Class XiEqModel Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 class XiEqModel
59 {
60 
61 protected:
62 
63  // Protected data
64 
65  //- Dictionary
67 
68  //- Thermo
70 
71  //- Turbulence
73 
74  //- Laminar burning velocity
75  const volScalarField& Su_;
76 
77 
78 public:
79 
80  //- Runtime type information
81  TypeName("XiEqModel");
82 
83 
84  // Declare run-time constructor selection table
85 
87  (
88  autoPtr,
89  XiEqModel,
90  dictionary,
91  (
92  const dictionary& XiEqProperties,
95  const volScalarField& Su
96  ),
97  (
98  XiEqProperties,
99  thermo,
100  turbulence,
101  Su
102  )
103  );
104 
105 
106  // Constructors
107 
108  //- Construct from components
109  XiEqModel
110  (
111  const dictionary& XiEqProperties,
112  const psiuReactionThermo& thermo,
113  const compressible::RASModel& turbulence,
114  const volScalarField& Su
115  );
116 
117  //- Disallow default bitwise copy construction
118  XiEqModel(const XiEqModel&);
119 
120 
121  // Selectors
122 
123  //- Return a reference to the selected XiEq model
124  static autoPtr<XiEqModel> New
125  (
126  const dictionary& XiEqProperties,
127  const psiuReactionThermo& thermo,
128  const compressible::RASModel& turbulence,
129  const volScalarField& Su
130  );
131 
132 
133  //- Destructor
134  virtual ~XiEqModel();
135 
136 
137  // Member Functions
138 
139  //- Return the flame-wrinkling XiEq
140  virtual tmp<volScalarField> XiEq() const
141  {
142  return turbulence_.muEff();
143  }
144 
145  //- Return the sub-grid Schelkin effect
146  tmp<volScalarField> calculateSchelkinEffect(const scalar) const;
147 
148  //- Update properties from given dictionary
149  virtual bool read(const dictionary& XiEqProperties) = 0;
150 
151  //- Write fields
152  void writeFields() const;
153 
154 
155  // Member Operators
156 
157  //- Disallow default bitwise assignment
158  void operator=(const XiEqModel&) = delete;
159 };
160 
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 } // End namespace Foam
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 #endif
169 
170 // ************************************************************************* //
dictionary XiEqModelCoeffs_
Dictionary.
Definition: XiEqModel.H:65
zeroField Su
Definition: alphaSuSp.H:1
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
tmp< volScalarField > calculateSchelkinEffect(const scalar) const
Return the sub-grid Schelkin effect.
declareRunTimeSelectionTable(autoPtr, XiEqModel, dictionary,(const dictionary &XiEqProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su),(XiEqProperties, thermo, turbulence, Su))
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
rhoReactionThermo & thermo
Definition: createFields.H:28
virtual ~XiEqModel()
Destructor.
XiEqModel(const dictionary &XiEqProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su)
Construct from components.
Base-class for all XiEq models used by the b-XiEq combustion model. The available models are : basicX...
Definition: XiEqModel.H:57
void operator=(const XiEqModel &)=delete
Disallow default bitwise assignment.
TypeName("XiEqModel")
Runtime type information.
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinkling XiEq.
Definition: XiEqModel.H:139
const psiuReactionThermo & thermo_
Thermo.
Definition: XiEqModel.H:68
Foam::psiuReactionThermo.
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
const volScalarField & Su_
Laminar burning velocity.
Definition: XiEqModel.H:74
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
static autoPtr< XiEqModel > New(const dictionary &XiEqProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su)
Return a reference to the selected XiEq model.
void writeFields() const
Write fields.
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
Namespace for OpenFOAM.
const compressible::RASModel & turbulence_
Turbulence.
Definition: XiEqModel.H:71