Gulders.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 "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 
50 (
51  const dictionary& dict,
52  const psiuReactionThermo& ct
53 )
54 :
56 
57  coeffsDict_(dict.subDict(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 
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 
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 
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 
158 (
159  const volScalarField& p,
160  const volScalarField& Tu,
161  const volScalarField& phi
162 ) const
163 {
165  (
166  new volScalarField
167  (
168  IOobject
169  (
170  "Su0",
171  p.time().timeName(),
172  p.db(),
175  false
176  ),
177  p.mesh(),
178  dimensionedScalar("Su0", dimVelocity, 0.0)
179  )
180  );
181 
182  volScalarField& Su0 = tSu0();
183 
184  forAll(Su0, celli)
185  {
186  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli], 0.0);
187  }
188 
189  forAll(Su0.boundaryField(), patchi)
190  {
191  forAll(Su0.boundaryField()[patchi], facei)
192  {
193  Su0.boundaryField()[patchi][facei] =
194  Su0pTphi
195  (
196  p.boundaryField()[patchi][facei],
197  Tu.boundaryField()[patchi][facei],
198  phi.boundaryField()[patchi][facei],
199  0.0
200  );
201  }
202  }
203 
204  return tSu0;
205 }
206 
207 
210 {
211  if (psiuReactionThermo_.composition().contains("ft"))
212  {
213  const volScalarField& ft = psiuReactionThermo_.composition().Y("ft");
214 
215  return Su0pTphi
216  (
217  psiuReactionThermo_.p(),
218  psiuReactionThermo_.Tu(),
220  (
221  psiuReactionThermo_.lookup("stoichiometricAirFuelMassRatio")
222  )*ft/max(1 - ft, SMALL)
223  );
224  }
225  else
226  {
227  return Su0pTphi
228  (
229  psiuReactionThermo_.p(),
230  psiuReactionThermo_.Tu(),
231  equivalenceRatio_
232  );
233  }
234 }
235 
236 
237 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::laminarFlameSpeedModels::Gulders::Gulders
Gulders(const Gulders &)
Construct as copy (not implemented)
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::laminarFlameSpeedModels::Gulders::~Gulders
virtual ~Gulders()
Destructor.
Definition: Gulders.C:69
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
scalar Su0pTphi(scalar p, scalar Tu, scalar phi, scalar Yres) const
Definition: Gulders.C:92
Foam::laminarFlameSpeedModels::Gulders::operator()
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Definition: Gulders.C:209
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::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::laminarFlameSpeed
Abstract class for laminar flame speed.
Definition: laminarFlameSpeed.H:58
Gulders.H
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::laminarFlameSpeedModels::Gulders::SuRef
scalar SuRef(scalar phi) const
Definition: Gulders.C:76