GuldersEGR.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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 
48 Foam::laminarFlameSpeedModels::GuldersEGR::GuldersEGR
49 (
50  const dictionary& dict,
51  const psiuReactionThermo& ct
52 )
53 :
54  laminarFlameSpeed(dict, ct),
55 
56  coeffsDict_(dict.optionalSubDict(typeName + "Coeffs").subDict(fuel_)),
57  W_(coeffsDict_.lookup<scalar>("W")),
58  eta_(coeffsDict_.lookup<scalar>("eta")),
59  xi_(coeffsDict_.lookup<scalar>("xi")),
60  f_(coeffsDict_.lookup<scalar>("f")),
61  alpha_(coeffsDict_.lookup<scalar>("alpha")),
62  beta_(coeffsDict_.lookup<scalar>("beta"))
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
73 
74 inline Foam::scalar Foam::laminarFlameSpeedModels::GuldersEGR::SuRef
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 
90 inline Foam::scalar Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
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 
106 Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
107 (
108  const volScalarField& p,
109  const volScalarField& Tu,
110  scalar phi
111 ) const
112 {
114  (
116  (
117  "Su0",
118  p.mesh(),
120  )
121  );
122 
123  volScalarField& Su0 = tSu0.ref();
124 
125  forAll(Su0, celli)
126  {
127  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi, 0.0);
128  }
129 
130  volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
131 
132  forAll(Su0Bf, patchi)
133  {
134  forAll(Su0Bf[patchi], facei)
135  {
136  Su0Bf[patchi][facei] =
137  Su0pTphi
138  (
139  p.boundaryField()[patchi][facei],
140  Tu.boundaryField()[patchi][facei],
141  phi,
142  0.0
143  );
144  }
145  }
146 
147  return tSu0;
148 }
149 
150 
152 Foam::laminarFlameSpeedModels::GuldersEGR::Su0pTphi
153 (
154  const volScalarField& p,
155  const volScalarField& Tu,
156  const volScalarField& phi,
157  const volScalarField& egr
158 ) const
159 {
161  (
163  (
164  "Su0",
165  p.mesh(),
167  )
168  );
169 
170  volScalarField& Su0 = tSu0.ref();
171 
172  forAll(Su0, celli)
173  {
174  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli], egr[celli]);
175  }
176 
177  volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
178 
179  forAll(Su0Bf, patchi)
180  {
181  forAll(Su0Bf[patchi], facei)
182  {
183  Su0Bf[patchi][facei] =
184  Su0pTphi
185  (
186  p.boundaryField()[patchi][facei],
187  Tu.boundaryField()[patchi][facei],
188  phi.boundaryField()[patchi][facei],
189  egr.boundaryField()[patchi][facei]
190  );
191  }
192  }
193 
194  return tSu0;
195 }
196 
197 
200 {
201  if
202  (
203  psiuReactionThermo_.composition().contains("ft")
204  && psiuReactionThermo_.composition().contains("egr")
205  )
206  {
207  return Su0pTphi
208  (
209  psiuReactionThermo_.p(),
210  psiuReactionThermo_.Tu(),
212  (
213  psiuReactionThermo_.lookup("stoichiometricAirFuelMassRatio")
214  )/
215  (
216  scalar(1)/psiuReactionThermo_.composition().Y("ft")
217  - scalar(1)
218  ),
219  psiuReactionThermo_.composition().Y("egr")
220  );
221  }
222  else
223  {
224  return Su0pTphi
225  (
226  psiuReactionThermo_.p(),
227  psiuReactionThermo_.Tu(),
228  equivalenceRatio_
229  );
230  }
231 }
232 
233 
234 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Definition: GuldersEGR.C:199
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Macros for easy insertion into run-time selection tables.
phi
Definition: pEqn.H:104
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:934
dimensionedScalar exp(const dimensionedScalar &ds)
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1001
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
Foam::psiuReactionThermo.
const Mesh & mesh() const
Return mesh.
scalar pRef
Definition: createFields.H:19
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
Abstract class for laminar flame speed.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
const dimensionSet dimVelocity