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-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::XiEqModel
26 
27 Description
28  Base class for equilibrium flame wrinkling \c XiEq models
29 
30 SourceFiles
31  XiEqModel.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef XiEqModel_H
36 #define XiEqModel_H
37 
38 #include "IOdictionary.H"
41 #include "runTimeSelectionTables.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class XiEqModel Declaration
50 \*---------------------------------------------------------------------------*/
51 
52 class XiEqModel
53 {
54 
55 protected:
56 
57  // Protected data
58 
59  //- Thermo
61 
62  //- Thermo-physical transport
64 
65  //- Turbulence
67 
68  //- Laminar burning velocity
69  const volScalarField& Su_;
70 
71  //- Update coefficients from given dictionary
72  virtual bool readCoeffs(const dictionary& dict) = 0;
73 
74 
75 public:
76 
77  //- Runtime type information
78  TypeName("XiEqModel");
79 
80 
81  // Declare run-time constructor selection table
82 
84  (
85  autoPtr,
86  XiEqModel,
87  dictionary,
88  (
89  const dictionary& dict,
92  const volScalarField& Su
93  ),
94  (
95  dict,
96  thermo,
97  turbulence,
98  Su
99  )
100  );
101 
102 
103  // Constructors
104 
105  //- Construct from components
106  XiEqModel
107  (
110  const volScalarField& Su
111  );
112 
113  //- Disallow default bitwise copy construction
114  XiEqModel(const XiEqModel&) = delete;
115 
116 
117  // Selectors
118 
119  //- Return a reference to the selected XiEq model
120  static autoPtr<XiEqModel> New
121  (
122  const dictionary& XiProperties,
125  const volScalarField& Su
126  );
127 
128 
129  //- Destructor
130  virtual ~XiEqModel();
131 
132 
133  // Member Functions
134 
135  //- Return the flame-wrinkling XiEq
136  virtual tmp<volScalarField> XiEq() const = 0;
137 
138  //- Return the flame diffusivity
139  virtual tmp<volScalarField> Db() const
140  {
141  return volScalarField::New
142  (
143  "Db",
145  );
146  }
147 
148  //- Update properties from the given XiProperties dictionary
149  bool read(const dictionary& XiProperties);
150 
151 
152  // Member Operators
153 
154  //- Disallow default bitwise assignment
155  void operator=(const XiEqModel&) = delete;
156 };
157 
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
161 } // End namespace Foam
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 #endif
166 
167 // ************************************************************************* //
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 equilibrium flame wrinkling XiEq models.
Definition: XiEqModel.H:52
const volScalarField & Su_
Laminar burning velocity.
Definition: XiEqModel.H:68
TypeName("XiEqModel")
Runtime type information.
virtual ~XiEqModel()
Destructor.
Definition: XiEqModel.C:63
const psiuMulticomponentThermo & thermo_
Thermo.
Definition: XiEqModel.H:59
static autoPtr< XiEqModel > New(const dictionary &XiProperties, const psiuMulticomponentThermo &thermo, const fluidThermoThermophysicalTransportModel &turbulence, const volScalarField &Su)
Return a reference to the selected XiEq model.
Definition: XiEqModelNew.C:31
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
Definition: XiEqModel.H:138
const compressibleMomentumTransportModel & turbulence_
Turbulence.
Definition: XiEqModel.H:65
void operator=(const XiEqModel &)=delete
Disallow default bitwise assignment.
virtual tmp< volScalarField > XiEq() const =0
Return the flame-wrinkling XiEq.
bool read(const dictionary &XiProperties)
Update properties from the given XiProperties dictionary.
Definition: XiEqModel.C:69
const fluidThermoThermophysicalTransportModel & thermoTransport_
Thermo-physical transport.
Definition: XiEqModel.H:62
XiEqModel(const psiuMulticomponentThermo &thermo, const fluidThermoThermophysicalTransportModel &turbulence, const volScalarField &Su)
Construct from components.
Definition: XiEqModel.C:48
virtual bool readCoeffs(const dictionary &dict)=0
Update coefficients from given dictionary.
Definition: XiEqModel.C:39
declareRunTimeSelectionTable(autoPtr, XiEqModel, dictionary,(const dictionary &dict, const psiuMulticomponentThermo &thermo, const fluidThermoThermophysicalTransportModel &turbulence, const volScalarField &Su),(dict, thermo, turbulence, Su))
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.
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 > DEff(const volScalarField &Yi) const =0
Effective mass diffusion coefficient.
Base-class for combustion fluid thermodynamic properties based on compressibility.
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
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