Gulder.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2012 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::XiEqModels::Gulder
28
29Description
30 Simple Gulder model for XiEq based on Gulders correlation
31 with a linear correction function to give a plausible profile for XiEq.
32
33SourceFiles
34 Gulder.C
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef Gulder_H
39#define Gulder_H
40
41#include "XiEqModel.H"
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45namespace Foam
46{
47namespace XiEqModels
48{
49
50/*---------------------------------------------------------------------------*\
51 Class Gulder Declaration
52\*---------------------------------------------------------------------------*/
54class Gulder
55:
56 public XiEqModel
57{
58 // Private data
59
60 //- Model constant
61 scalar XiEqCoef_;
62
63 //- Minimum laminar burning velocity
64 const dimensionedScalar SuMin_;
65
66 //- Schelkin effect Model constant
67 scalar uPrimeCoef_;
68
69 //- Use sub-grid Schelkin effect
70 bool subGridSchelkin_;
71
72
73 // Private Member Functions
74
75 //- No copy construct
76 Gulder(const Gulder&) = delete;
77
78 //- No copy assignment
79 void operator=(const Gulder&) = delete;
80
81
82public:
83
84 //- Runtime type information
85 TypeName("Gulder");
86
87
88 // Constructors
89
90 //- Construct from components
91 Gulder
92 (
93 const dictionary& XiEqProperties,
95 const compressible::RASModel& turbulence,
96 const volScalarField& Su
97 );
98
99
100 //- Destructor
101 virtual ~Gulder();
102
103
104 // Member Functions
105
106 //- Return the flame-wrinkling XiEq
107 virtual tmp<volScalarField> XiEq() const;
108
109 //- Update properties from given dictionary
110 virtual bool read(const dictionary& XiEqProperties);
111
112};
113
114
115// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
116
117} // End namespace XiEqModels
118} // End namespace Foam
119
120// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121
122#endif
123
124// ************************************************************************* //
Base-class for all XiEq models used by the b-XiEq combustion model. The available models are : basicX...
Definition: XiEqModel.H:60
Simple Gulder model for XiEq based on Gulders correlation with a linear correction function to give a...
Definition: Gulder.H:56
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinkling XiEq.
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
Gulder(const dictionary &XiEqProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su)
Construct from components.
virtual ~Gulder()
Destructor.
TypeName("Gulder")
Runtime type information.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Foam::psiuReactionThermo.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
zeroField Su
Definition: alphaSuSp.H:1
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
Namespace for OpenFOAM.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73