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 | Copyright (C) 2011-2015 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 \*---------------------------------------------------------------------------*/
25 
26 #include "GuldersEGR.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace laminarFlameSpeedModels
34 {
35  defineTypeNameAndDebug(GuldersEGR, 0);
36 
38  (
39  laminarFlameSpeed,
40  GuldersEGR,
41  dictionary
42  );
43 }
44 }
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const dictionary& dict,
51  const psiuReactionThermo& ct
52 )
53 :
55 
56  coeffsDict_(dict.subDict(typeName + "Coeffs").subDict(fuel_)),
57  W_(readScalar(coeffsDict_.lookup("W"))),
58  eta_(readScalar(coeffsDict_.lookup("eta"))),
59  xi_(readScalar(coeffsDict_.lookup("xi"))),
60  f_(readScalar(coeffsDict_.lookup("f"))),
61  alpha_(readScalar(coeffsDict_.lookup("alpha"))),
62  beta_(readScalar(coeffsDict_.lookup("beta")))
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
73 
75 (
76  scalar phi
77 ) const
78 {
79  if (phi > SMALL)
80  {
81  return W_*pow(phi, eta_)*exp(-xi_*sqr(phi - 1.075));
82  }
83  else
84  {
85  return 0.0;
86  }
87 }
88 
89 
91 (
92  scalar p,
93  scalar Tu,
94  scalar phi,
95  scalar Yres
96 ) const
97 {
98  static const scalar Tref = 300.0;
99  static const scalar pRef = 1.013e5;
100 
101  return SuRef(phi)*pow((Tu/Tref), alpha_)*pow((p/pRef), beta_)*(1 - f_*Yres);
102 }
103 
104 
107 (
108  const volScalarField& p,
109  const volScalarField& Tu,
110  scalar phi
111 ) const
112 {
114  (
115  new volScalarField
116  (
117  IOobject
118  (
119  "Su0",
120  p.time().timeName(),
121  p.db(),
124  false
125  ),
126  p.mesh(),
127  dimensionedScalar("Su0", dimVelocity, 0.0)
128  )
129  );
130 
131  volScalarField& Su0 = tSu0();
132 
133  forAll(Su0, celli)
134  {
135  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi, 0.0);
136  }
137 
138  forAll(Su0.boundaryField(), patchi)
139  {
140  forAll(Su0.boundaryField()[patchi], facei)
141  {
142  Su0.boundaryField()[patchi][facei] =
143  Su0pTphi
144  (
145  p.boundaryField()[patchi][facei],
146  Tu.boundaryField()[patchi][facei],
147  phi,
148  0.0
149  );
150  }
151  }
152 
153  return tSu0;
154 }
155 
156 
159 (
160  const volScalarField& p,
161  const volScalarField& Tu,
162  const volScalarField& phi,
163  const volScalarField& egr
164 ) const
165 {
167  (
168  new volScalarField
169  (
170  IOobject
171  (
172  "Su0",
173  p.time().timeName(),
174  p.db(),
177  false
178  ),
179  p.mesh(),
180  dimensionedScalar("Su0", dimVelocity, 0.0)
181  )
182  );
183 
184  volScalarField& Su0 = tSu0();
185 
186  forAll(Su0, celli)
187  {
188  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli], egr[celli]);
189  }
190 
191  forAll(Su0.boundaryField(), patchi)
192  {
193  forAll(Su0.boundaryField()[patchi], facei)
194  {
195  Su0.boundaryField()[patchi][facei] =
196  Su0pTphi
197  (
198  p.boundaryField()[patchi][facei],
199  Tu.boundaryField()[patchi][facei],
200  phi.boundaryField()[patchi][facei],
201  egr.boundaryField()[patchi][facei]
202  );
203  }
204  }
205 
206  return tSu0;
207 }
208 
209 
212 {
213  if
214  (
215  psiuReactionThermo_.composition().contains("ft")
216  && psiuReactionThermo_.composition().contains("egr")
217  )
218  {
219  return Su0pTphi
220  (
221  psiuReactionThermo_.p(),
222  psiuReactionThermo_.Tu(),
224  (
225  psiuReactionThermo_.lookup("stoichiometricAirFuelMassRatio")
226  )/
227  (
228  scalar(1)/psiuReactionThermo_.composition().Y("ft")
229  - scalar(1)
230  ),
231  psiuReactionThermo_.composition().Y("egr")
232  );
233  }
234  else
235  {
236  return Su0pTphi
237  (
238  psiuReactionThermo_.p(),
239  psiuReactionThermo_.Tu(),
240  equivalenceRatio_
241  );
242  }
243 }
244 
245 
246 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
p
p
Definition: pEqn.H:62
Foam::laminarFlameSpeedModels::defineTypeNameAndDebug
defineTypeNameAndDebug(constant, 0)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::psiuReactionThermo
Foam::psiuReactionThermo.
Definition: psiuReactionThermo.H:49
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
scalar Su0pTphi(scalar p, scalar Tu, scalar phi, scalar Yres) const
Definition: GuldersEGR.C:91
Foam::laminarFlameSpeedModels::GuldersEGR::GuldersEGR
GuldersEGR(const GuldersEGR &)
Construct as copy (not implemented)
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
GuldersEGR.H
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
pRef
scalar pRef
Definition: createFields.H:19
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::laminarFlameSpeedModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
readScalar
#define readScalar
Definition: doubleScalar.C:38
Foam::laminarFlameSpeedModels::GuldersEGR::SuRef
scalar SuRef(scalar phi) const
Definition: GuldersEGR.C:75
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::laminarFlameSpeed
Abstract class for laminar flame speed.
Definition: laminarFlameSpeed.H:58
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::laminarFlameSpeedModels::GuldersEGR::~GuldersEGR
virtual ~GuldersEGR()
Destructor.
Definition: GuldersEGR.C:68
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::laminarFlameSpeedModels::GuldersEGR::operator()
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Definition: GuldersEGR.C:211