Gulder.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-2023 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::XiEqModels::Gulder
26 
27 Description
28  Simple Gulder model for XiEq based on Gulders correlation
29  with a linear correction function to give a plausible profile for XiEq.
30 
31 SourceFiles
32  Gulder.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef Gulder_H
37 #define Gulder_H
38 
39 #include "XiEqModel.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace XiEqModels
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class Gulder Declaration
50 \*---------------------------------------------------------------------------*/
51 
52 class Gulder
53 :
54  public XiEqModel
55 {
56  // Private Data
57 
58  //- Model constant
59  scalar XiEqCoef_;
60 
61  //- Minimum laminar burning velocity
62  const dimensionedScalar SuMin_;
63 
64  //- Schelkin effect Model constant
65  scalar uPrimeCoef_;
66 
67  //- Use sub-grid Schelkin effect
68  bool subGridSchelkin_;
69 
70 
71 public:
72 
73  //- Runtime type information
74  TypeName("Gulder");
75 
76 
77  // Constructors
78 
79  //- Construct from components
80  Gulder
81  (
82  const dictionary& XiEqProperties,
85  const volScalarField& Su
86  );
87 
88  //- Disallow default bitwise copy construction
89  Gulder(const Gulder&) = delete;
90 
91 
92  //- Destructor
93  virtual ~Gulder();
94 
95 
96  // Member Functions
97 
98  //- Return the flame-wrinkling XiEq
99  virtual tmp<volScalarField> XiEq() const;
100 
101  //- Update properties from given dictionary
102  virtual bool read(const dictionary& XiEqProperties);
103 
104 
105  // Member Operators
106 
107  //- Disallow default bitwise assignment
108  void operator=(const Gulder&) = delete;
109 };
110 
111 
112 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
113 
114 } // End namespace XiEqModels
115 } // End namespace Foam
116 
117 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
118 
119 #endif
120 
121 // ************************************************************************* //
Generic GeometricField class.
Base-class for all XiEq models used by the b-XiEq combustion model. The available models are : basicX...
Definition: XiEqModel.H:58
Simple Gulder model for XiEq based on Gulders correlation with a linear correction function to give a...
Definition: Gulder.H:54
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
Definition: Gulder.C:97
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinkling XiEq.
Definition: Gulder.C:70
void operator=(const Gulder &)=delete
Disallow default bitwise assignment.
Gulder(const dictionary &XiEqProperties, const psiuMulticomponentThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su)
Construct from components.
Definition: Gulder.C:44
TypeName("Gulder")
Runtime type information.
virtual ~Gulder()
Destructor.
Definition: Gulder.C:64
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Base-class for combustion fluid thermodynamic properties based on compressibility.
A class for managing temporary objects.
Definition: tmp.H:55
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
Namespace for OpenFOAM.
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
fluidMulticomponentThermo & thermo
Definition: createFields.H:31