GuldersEGR.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 "GuldersEGR.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace laminarFlameSpeedModels
36{
38
40 (
44 );
45}
46}
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
51(
52 const dictionary& dict,
53 const psiuReactionThermo& ct
54)
55:
57
58 coeffsDict_(dict.optionalSubDict(typeName + "Coeffs").subDict(fuel_)),
59 W_(coeffsDict_.get<scalar>("W")),
60 eta_(coeffsDict_.get<scalar>("eta")),
61 xi_(coeffsDict_.get<scalar>("xi")),
62 f_(coeffsDict_.get<scalar>("f")),
63 alpha_(coeffsDict_.get<scalar>("alpha")),
64 beta_(coeffsDict_.get<scalar>("beta"))
65{}
66
67
68// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69
71{}
72
73
74// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
75
76inline Foam::scalar Foam::laminarFlameSpeedModels::GuldersEGR::SuRef
77(
78 scalar phi
79) const
80{
81 if (phi > SMALL)
82 {
83 return W_*pow(phi, eta_)*exp(-xi_*sqr(phi - 1.075));
84 }
85 else
86 {
87 return 0.0;
88 }
89}
90
91
92inline Foam::scalar Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
93(
94 scalar p,
95 scalar Tu,
96 scalar phi,
97 scalar Yres
98) const
99{
100 static const scalar Tref = 300.0;
101 static const scalar pRef = 1.013e5;
102
103 return SuRef(phi)*pow((Tu/Tref), alpha_)*pow((p/pRef), beta_)*(1 - f_*Yres);
104}
105
106
108Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
109(
110 const volScalarField& p,
111 const volScalarField& Tu,
112 scalar phi
113) const
114{
115 tmp<volScalarField> tSu0
116 (
118 (
119 IOobject
120 (
121 "Su0",
122 p.time().timeName(),
123 p.db(),
126 false
127 ),
128 p.mesh(),
130 )
131 );
132
133 volScalarField& Su0 = tSu0.ref();
134
135 forAll(Su0, celli)
136 {
137 Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi, 0.0);
138 }
139
140 volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
141
142 forAll(Su0Bf, patchi)
143 {
144 forAll(Su0Bf[patchi], facei)
145 {
146 Su0Bf[patchi][facei] =
147 Su0pTphi
148 (
149 p.boundaryField()[patchi][facei],
150 Tu.boundaryField()[patchi][facei],
151 phi,
152 0.0
153 );
154 }
155 }
156
157 return tSu0;
158}
159
160
162Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
163(
164 const volScalarField& p,
165 const volScalarField& Tu,
166 const volScalarField& phi,
167 const volScalarField& egr
168) const
169{
170 tmp<volScalarField> tSu0
171 (
173 (
174 IOobject
175 (
176 "Su0",
177 p.time().timeName(),
178 p.db(),
181 false
182 ),
183 p.mesh(),
185 )
186 );
187
188 volScalarField& Su0 = tSu0.ref();
189
190 forAll(Su0, celli)
191 {
192 Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli], egr[celli]);
193 }
194
195 volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
196
197 forAll(Su0Bf, patchi)
198 {
199 forAll(Su0Bf[patchi], facei)
200 {
201 Su0Bf[patchi][facei] =
202 Su0pTphi
203 (
204 p.boundaryField()[patchi][facei],
205 Tu.boundaryField()[patchi][facei],
206 phi.boundaryField()[patchi][facei],
207 egr.boundaryField()[patchi][facei]
208 );
209 }
210 }
211
212 return tSu0;
213}
214
215
218{
219 if
220 (
221 psiuReactionThermo_.composition().contains("ft")
222 && psiuReactionThermo_.composition().contains("egr")
223 )
224 {
225 return Su0pTphi
226 (
227 psiuReactionThermo_.p(),
228 psiuReactionThermo_.Tu(),
230 (
231 "stoichiometricAirFuelMassRatio", dimless, psiuReactionThermo_
232 )/
233 (
234 scalar(1)/psiuReactionThermo_.composition().Y("ft")
235 - scalar(1)
236 ),
237 psiuReactionThermo_.composition().Y("egr")
238 );
239 }
240 else
241 {
242 return Su0pTphi
243 (
244 psiuReactionThermo_.p(),
245 psiuReactionThermo_.Tu(),
246 equivalenceRatio_
247 );
248 }
249}
250
251
252// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Laminar flame speed obtained from Gulder's correlation with EGR modelling.
Definition: GuldersEGR.H:55
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Definition: GuldersEGR.C:217
Abstract class for laminar flame speed.
Foam::psiuReactionThermo.
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimVelocity
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333