OpenFOAM
4.1
The OpenFOAM Foundation
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Modules
Pages
applications
solvers
combustion
PDRFoam
PDRModels
XiEqModels
basicXiSubXiEq
basicXiSubXiEq.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) 2011-2012 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::basicSubGrid
26
27
Description
28
Basic sub-grid obstacle flame-wrinking enhancement factor model.
29
Details supplied by J Puttock 2/7/06.
30
31
<b> Sub-grid flame area generation </b>
32
33
\f$ n = N - \hat{\dwea{\vec{U}}}.n_{s}.\hat{\dwea{\vec{U}}} \f$
34
\f$ n_{r} = \sqrt{n} \f$
35
36
where:
37
38
\f$ \hat{\dwea{\vec{U}}} = \dwea{\vec{U}} / \vert \dwea{\vec{U}}
39
\vert \f$
40
41
\f$ b = \hat{\dwea{\vec{U}}}.B.\hat{\dwea{\vec{U}}} / n_{r} \f$
42
43
where:
44
45
\f$ B \f$ is the file "B".
46
47
\f$ N \f$ is the file "N".
48
49
\f$ n_{s} \f$ is the file "ns".
50
51
The flame area enhancement factor \f$ \Xi_{sub} \f$ is expected to
52
approach:
53
54
\f[
55
\Xi_{{sub}_{eq}} =
56
1 + max(2.2 \sqrt{b}, min(0.34 \frac{\vert \dwea{\vec{U}}
57
\vert}{{\vec{U}}^{'}}, 1.6)) \times min(\frac{n}{4}, 1)
58
\f]
59
60
61
SourceFiles
62
basicSubGrid.C
63
64
\*---------------------------------------------------------------------------*/
65
66
#ifndef basicSubGrid_H
67
#define basicSubGrid_H
68
69
#include "
XiEqModel.H
"
70
71
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72
73
namespace
Foam
74
{
75
namespace
XiEqModels
76
{
77
78
/*---------------------------------------------------------------------------*\
79
Class basicSubGrid Declaration
80
\*---------------------------------------------------------------------------*/
81
82
class
basicSubGrid
83
:
84
public
XiEqModel
85
{
86
// Private data
87
88
//- tblock
89
volSymmTensorField
B_;
90
91
//- Equilibrium Xi model due to turbulence
92
autoPtr<XiEqModel>
XiEqModel_;
93
94
95
// Private Member Functions
96
97
//- Disallow copy construct
98
basicSubGrid
(
const
basicSubGrid
&);
99
100
//- Disallow default bitwise assignment
101
void
operator=(
const
basicSubGrid
&);
102
103
104
public
:
105
106
//- Runtime type information
107
TypeName
(
"basicSubGrid"
);
108
109
110
// Constructors
111
112
//- Construct from components
113
basicSubGrid
114
(
115
const
dictionary
& XiEqProperties,
116
const
psiuReactionThermo
&
thermo
,
117
const
compressible::RASModel
&
turbulence
,
118
const
volScalarField
&
Su
119
);
120
121
122
//- Destructor
123
virtual
~basicSubGrid
();
124
125
126
// Member Functions
127
128
//- Return the flame-wrinking XiEq
129
virtual
tmp<volScalarField>
XiEq
()
const
;
130
131
//- Update properties from given dictionary
132
virtual
bool
read
(
const
dictionary
& XiEqProperties);
133
};
134
135
136
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137
138
}
// End namespace XiEqModels
139
}
// End namespace Foam
140
141
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142
143
#endif
144
145
// ************************************************************************* //
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition:
createFields.H:23
Foam::XiEqModels::basicSubGrid::TypeName
TypeName("basicSubGrid")
Runtime type information.
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition:
dictionary.H:137
Foam::XiEqModels::basicSubGrid
Basic sub-grid obstacle flame-wrinking enhancement factor model. Details supplied by J Puttock 2/7/06...
Definition:
basicXiSubXiEq.H:81
Foam::GeometricField< symmTensor, fvPatchField, volMesh >
XiEqModel.H
Foam::XiEqModel
Base-class for all XiEq models used by the b-XiEq combustion model. The available models are : basicX...
Definition:
XiEqModel.H:57
Foam::XiEqModels::basicSubGrid::XiEq
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinking XiEq.
Foam::XiEqModels::basicSubGrid::~basicSubGrid
virtual ~basicSubGrid()
Destructor.
Foam::psiuReactionThermo
Foam::psiuReactionThermo.
Definition:
psiuReactionThermo.H:49
Foam::thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::XiEqModels::basicSubGrid::read
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
Foam::compressible::RASModel
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
Definition:
turbulentFluidThermoModel.H:62
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition:
PtrList.H:53
Foam::tmp
A class for managing temporary objects.
Definition:
PtrList.H:54
Foam::fvc::Su
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition:
fvcSup.C:44
Foam
Namespace for OpenFOAM.
Definition:
combustionModel.C:30
Generated by
1.8.11