Gulders.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-2018 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 "Gulders.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace laminarFlameSpeedModels
34 {
35  defineTypeNameAndDebug(Gulders, 0);
36 
38  (
39  laminarFlameSpeed,
40  Gulders,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 Foam::laminarFlameSpeedModels::Gulders::Gulders
50 (
51  const dictionary& dict,
52  const psiuReactionThermo& ct
53 )
54 :
55  laminarFlameSpeed(dict, ct),
56 
57  coeffsDict_(dict.optionalSubDict(typeName + "Coeffs").subDict(fuel_)),
58  W_(readScalar(coeffsDict_.lookup("W"))),
59  eta_(readScalar(coeffsDict_.lookup("eta"))),
60  xi_(readScalar(coeffsDict_.lookup("xi"))),
61  f_(readScalar(coeffsDict_.lookup("f"))),
62  alpha_(readScalar(coeffsDict_.lookup("alpha"))),
63  beta_(readScalar(coeffsDict_.lookup("beta")))
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
68 
70 {}
71 
72 
73 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
74 
75 inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::SuRef
76 (
77  scalar phi
78 ) const
79 {
80  if (phi > small)
81  {
82  return W_*pow(phi, eta_)*exp(-xi_*sqr(phi - 1.075));
83  }
84  else
85  {
86  return 0.0;
87  }
88 }
89 
90 
91 inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
92 (
93  scalar p,
94  scalar Tu,
95  scalar phi,
96  scalar Yres
97 ) const
98 {
99  static const scalar Tref = 300.0;
100  static const scalar pRef = 1.013e5;
101 
102  return SuRef(phi)*pow((Tu/Tref), alpha_)*pow((p/pRef), beta_)*(1 - f_*Yres);
103 }
104 
105 
106 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
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.ref();
132 
133  forAll(Su0, celli)
134  {
135  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi, 0.0);
136  }
137 
138  volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
139 
140  forAll(Su0Bf, patchi)
141  {
142  forAll(Su0Bf[patchi], facei)
143  {
144  Su0Bf[patchi][facei] =
145  Su0pTphi
146  (
147  p.boundaryField()[patchi][facei],
148  Tu.boundaryField()[patchi][facei],
149  phi,
150  0.0
151  );
152  }
153  }
154 
155  return tSu0;
156 }
157 
158 
159 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
160 (
161  const volScalarField& p,
162  const volScalarField& Tu,
163  const volScalarField& phi
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.ref();
185 
186  forAll(Su0, celli)
187  {
188  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli], 0.0);
189  }
190 
191  volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
192 
193  forAll(Su0Bf, patchi)
194  {
195  forAll(Su0Bf[patchi], facei)
196  {
197  Su0Bf[patchi][facei] =
198  Su0pTphi
199  (
200  p.boundaryField()[patchi][facei],
201  Tu.boundaryField()[patchi][facei],
202  phi.boundaryField()[patchi][facei],
203  0.0
204  );
205  }
206  }
207 
208  return tSu0;
209 }
210 
211 
214 {
215  if (psiuReactionThermo_.composition().contains("ft"))
216  {
217  const volScalarField& ft = psiuReactionThermo_.composition().Y("ft");
218 
219  return Su0pTphi
220  (
221  psiuReactionThermo_.p(),
222  psiuReactionThermo_.Tu(),
224  (
225  psiuReactionThermo_.lookup("stoichiometricAirFuelMassRatio")
226  )*ft/max(1 - ft, small)
227  );
228  }
229  else
230  {
231  return Su0pTphi
232  (
233  psiuReactionThermo_.p(),
234  psiuReactionThermo_.Tu(),
235  equivalenceRatio_
236  );
237  }
238 }
239 
240 
241 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
surfaceScalarField & phi
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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 word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
Macros for easy insertion into run-time selection tables.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
dimensionedScalar exp(const dimensionedScalar &ds)
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:759
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Definition: Gulders.C:213
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
Foam::psiuReactionThermo.
virtual ~Gulders()
Destructor.
Definition: Gulders.C:69
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.
const Time & time() const
Return time.
Definition: IOobject.C:367
A class for managing temporary objects.
Definition: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:361
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.
const dimensionSet dimVelocity