Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
applications
solvers
combustion
PDRFoam
XiModels
algebraic
algebraic.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 | Copyright (C) 1991-2010 OpenCFD Ltd.
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::XiModels::algebraic
26
27
Description
28
Simple algebraic model for Xi based on Gulders correlation
29
with a linear correction function to give a plausible profile for Xi.
30
See report TR/HGW/10 for details on the Weller two equations model.
31
See \link XiModel.H \endlink for more details on flame wrinkling modelling.
32
33
SourceFiles
34
algebraic.C
35
36
\*---------------------------------------------------------------------------*/
37
38
#ifndef algebraic_H
39
#define algebraic_H
40
41
#include "
../XiModel/XiModel.H
"
42
#include "
../XiEqModels/XiEqModel/XiEqModel.H
"
43
#include "
../XiGModels/XiGModel/XiGModel.H
"
44
45
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47
namespace
Foam
48
{
49
namespace
XiModels
50
{
51
52
/*---------------------------------------------------------------------------*\
53
Class algebraic Declaration
54
\*---------------------------------------------------------------------------*/
55
56
class
algebraic
57
:
58
public
XiModel
59
{
60
// Private data
61
62
scalar XiShapeCoef;
63
64
autoPtr<XiEqModel>
XiEqModel_;
65
autoPtr<XiGModel>
XiGModel_;
66
67
68
// Private Member Functions
69
70
//- Disallow copy construct
71
algebraic
(
const
algebraic
&);
72
73
//- Disallow default bitwise assignment
74
void
operator=(
const
algebraic
&);
75
76
77
public
:
78
79
//- Runtime type information
80
TypeName
(
"algebraic"
);
81
82
83
// Constructors
84
85
//- Construct from components
86
algebraic
87
(
88
const
dictionary
& XiProperties,
89
const
hhuCombustionThermo
&
thermo
,
90
const
compressible::RASModel
&
turbulence
,
91
const
volScalarField
&
Su
,
92
const
volScalarField
&
rho
,
93
const
volScalarField
&
b
,
94
const
surfaceScalarField
&
phi
95
);
96
97
98
// Destructor
99
100
virtual
~algebraic
();
101
102
103
// Member Functions
104
105
//- Return the flame diffusivity
106
virtual
tmp<volScalarField>
Db
()
const
;
107
108
//- Correct the flame-wrinking Xi
109
virtual
void
correct
();
110
111
//- Update properties from given dictionary
112
virtual
bool
read
(
const
dictionary
& XiProperties);
113
};
114
115
116
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
117
118
}
// End namespace XiModels
119
}
// End namespace Foam
120
121
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122
123
#endif
124
125
// ************************ vim: set sw=4 sts=4 et: ************************ //