Loading...
Searching...
No Matches
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"
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33namespace Foam
34{
36{
38
40 (
44 );
45}
46}
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
50Foam::laminarFlameSpeedModels::GuldersEGR::GuldersEGR
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
107Foam::tmp<Foam::volScalarField>
108Foam::laminarFlameSpeedModels::GuldersEGR::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>
151Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
152(
153 const volScalarField& p,
154 const volScalarField& Tu,
155 const volScalarField& phi,
156 const volScalarField& egr
157) const
158{
159 auto tSu0 = volScalarField::New
160 (
161 "Su0",
163 p.mesh(),
165 );
166 auto& Su0 = tSu0.ref();
167
168 forAll(Su0, celli)
169 {
170 Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli], egr[celli]);
171 }
172
173 volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
174
175 forAll(Su0Bf, patchi)
176 {
177 forAll(Su0Bf[patchi], facei)
178 {
179 Su0Bf[patchi][facei] =
180 Su0pTphi
181 (
182 p.boundaryField()[patchi][facei],
183 Tu.boundaryField()[patchi][facei],
184 phi.boundaryField()[patchi][facei],
185 egr.boundaryField()[patchi][facei]
186 );
187 }
189
190 return tSu0;
191}
192
193
196{
197 if
198 (
199 psiuReactionThermo_.composition().contains("ft")
200 && psiuReactionThermo_.composition().contains("egr")
201 )
202 {
203 return Su0pTphi
204 (
205 psiuReactionThermo_.p(),
206 psiuReactionThermo_.Tu(),
208 (
209 "stoichiometricAirFuelMassRatio", dimless, psiuReactionThermo_
210 )/
211 (
212 scalar(1)/psiuReactionThermo_.composition().Y("ft")
213 - scalar(1)
214 ),
215 psiuReactionThermo_.composition().Y("egr")
216 );
217 }
218 else
219 {
220 return Su0pTphi
221 (
222 psiuReactionThermo_.p(),
223 psiuReactionThermo_.Tu(),
224 equivalenceRatio_
225 );
226 }
227}
228
229
230// ************************************************************************* //
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 with EGR modelling.
Definition GuldersEGR.H:51
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Definition GuldersEGR.C:188
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.
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