Loading...
Searching...
No Matches
Gulders.C
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-2017 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
26\*---------------------------------------------------------------------------*/
27
28#include "Gulders.H"
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33namespace Foam
34{
36{
38
40 (
42 Gulders,
44 );
45}
46}
47
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
51Foam::laminarFlameSpeedModels::Gulders::Gulders
52(
53 const dictionary& dict,
54 const psiuReactionThermo& ct
55)
56:
58
59 coeffsDict_(dict.optionalSubDict(typeName + "Coeffs").subDict(fuel_)),
60 W_(coeffsDict_.get<scalar>("W")),
61 eta_(coeffsDict_.get<scalar>("eta")),
62 xi_(coeffsDict_.get<scalar>("xi")),
63 f_(coeffsDict_.get<scalar>("f")),
64 alpha_(coeffsDict_.get<scalar>("alpha")),
65 beta_(coeffsDict_.get<scalar>("beta"))
66{}
67
68
69// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
70
72{}
73
74
75// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
76
77inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::SuRef
78(
79 scalar phi
80) const
81{
82 if (phi > SMALL)
83 {
84 return W_*pow(phi, eta_)*exp(-xi_*sqr(phi - 1.075));
85 }
86 else
87 {
88 return 0.0;
89 }
90}
91
92
93inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
94(
95 scalar p,
96 scalar Tu,
97 scalar phi,
98 scalar Yres
99) const
100{
101 static const scalar Tref = 300.0;
102 static const scalar pRef = 1.013e5;
103
104 return SuRef(phi)*pow((Tu/Tref), alpha_)*pow((p/pRef), beta_)*(1 - f_*Yres);
105}
106
107
108Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
109(
110 const volScalarField& p,
111 const volScalarField& Tu,
112 scalar phi
113) const
114{
115 auto tSu0 = volScalarField::New
116 (
117 "Su0",
119 p.mesh(),
121 );
122 auto& Su0 = tSu0.ref();
123
124 forAll(Su0, celli)
125 {
126 Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi, 0.0);
127 }
128
129 volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
130
131 forAll(Su0Bf, patchi)
132 {
133 forAll(Su0Bf[patchi], facei)
134 {
135 Su0Bf[patchi][facei] =
136 Su0pTphi
137 (
138 p.boundaryField()[patchi][facei],
139 Tu.boundaryField()[patchi][facei],
140 phi,
141 0.0
142 );
143 }
144 }
145
146 return tSu0;
147}
148
149
150Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
151(
152 const volScalarField& p,
153 const volScalarField& Tu,
154 const volScalarField& phi
155) const
156{
157 auto tSu0 = volScalarField::New
158 (
159 "Su0",
161 p.mesh(),
163 );
164 auto& Su0 = tSu0.ref();
165
166 forAll(Su0, celli)
167 {
168 Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli], 0.0);
169 }
170
171 volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
172
173 forAll(Su0Bf, patchi)
174 {
175 forAll(Su0Bf[patchi], facei)
176 {
177 Su0Bf[patchi][facei] =
178 Su0pTphi
179 (
180 p.boundaryField()[patchi][facei],
181 Tu.boundaryField()[patchi][facei],
182 phi.boundaryField()[patchi][facei],
183 0.0
184 );
185 }
187
188 return tSu0;
189}
190
191
194{
195 if (psiuReactionThermo_.composition().contains("ft"))
196 {
197 const volScalarField& ft = psiuReactionThermo_.composition().Y("ft");
198
199 return Su0pTphi
200 (
201 psiuReactionThermo_.p(),
202 psiuReactionThermo_.Tu(),
204 (
205 "stoichiometricAirFuelMassRatio", dimless, psiuReactionThermo_
206 )*ft/max(1 - ft, SMALL)
207 );
208 }
209 else
210 {
211 return Su0pTphi
212 (
213 psiuReactionThermo_.p(),
214 psiuReactionThermo_.Tu(),
215 equivalenceRatio_
216 );
217 }
218}
219
220
221// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
@ NO_REGISTER
Do not request registration (bool: false).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Laminar flame speed obtained from Gulder's correlation.
Definition Gulders.H:51
virtual ~Gulders()
Destructor.
Definition Gulders.C:64
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Definition Gulders.C:186
Abstract class for laminar flame speed.
const psiuReactionThermo & psiuReactionThermo_
scalar equivalenceRatio_
Equivalence ratio of a homogeneous mixture.
Foam::psiuReactionThermo.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
Namespace for laminar flame speed models.
Definition constant.C:29
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299