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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 "Gulders.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace laminarFlameSpeedModels
36 {
37  defineTypeNameAndDebug(Gulders, 0);
38 
40  (
41  laminarFlameSpeed,
42  Gulders,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 Foam::laminarFlameSpeedModels::Gulders::Gulders
52 (
53  const dictionary& dict,
54  const psiuReactionThermo& ct
55 )
56 :
58 
59  coeffsDict_(dict.optionalSubDict(typeName + "Coeffs").subDict(fuel_)),
60  W_(coeffsDict_.get<scalar>("W")),
61  eta_(coeffsDict_.get<scalar>("eta")),
62  xi_(coeffsDict_.get<scalar>("xi")),
63  f_(coeffsDict_.get<scalar>("f")),
64  alpha_(coeffsDict_.get<scalar>("alpha")),
65  beta_(coeffsDict_.get<scalar>("beta"))
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
70 
72 {}
73 
74 
75 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
76 
77 inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::SuRef
78 (
79  scalar phi
80 ) const
81 {
82  if (phi > SMALL)
83  {
84  return W_*pow(phi, eta_)*exp(-xi_*sqr(phi - 1.075));
85  }
86  else
87  {
88  return 0.0;
89  }
90 }
91 
92 
93 inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
94 (
95  scalar p,
96  scalar Tu,
97  scalar phi,
98  scalar Yres
99 ) const
100 {
101  static const scalar Tref = 300.0;
102  static const scalar pRef = 1.013e5;
103 
104  return SuRef(phi)*pow((Tu/Tref), alpha_)*pow((p/pRef), beta_)*(1 - f_*Yres);
105 }
106 
107 
108 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
109 (
110  const volScalarField& p,
111  const volScalarField& Tu,
112  scalar phi
113 ) const
114 {
115  tmp<volScalarField> tSu0
116  (
117  new volScalarField
118  (
119  IOobject
120  (
121  "Su0",
122  p.time().timeName(),
123  p.db(),
126  false
127  ),
128  p.mesh(),
130  )
131  );
132 
133  volScalarField& Su0 = tSu0.ref();
134 
135  forAll(Su0, celli)
136  {
137  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi, 0.0);
138  }
139 
140  volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
141 
142  forAll(Su0Bf, patchi)
143  {
144  forAll(Su0Bf[patchi], facei)
145  {
146  Su0Bf[patchi][facei] =
147  Su0pTphi
148  (
149  p.boundaryField()[patchi][facei],
150  Tu.boundaryField()[patchi][facei],
151  phi,
152  0.0
153  );
154  }
155  }
156 
157  return tSu0;
158 }
159 
160 
161 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
162 (
163  const volScalarField& p,
164  const volScalarField& Tu,
165  const volScalarField& phi
166 ) const
167 {
168  tmp<volScalarField> tSu0
169  (
170  new volScalarField
171  (
172  IOobject
173  (
174  "Su0",
175  p.time().timeName(),
176  p.db(),
179  false
180  ),
181  p.mesh(),
183  )
184  );
185 
186  volScalarField& Su0 = tSu0.ref();
187 
188  forAll(Su0, celli)
189  {
190  Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli], 0.0);
191  }
192 
193  volScalarField::Boundary& Su0Bf = Su0.boundaryFieldRef();
194 
195  forAll(Su0Bf, patchi)
196  {
197  forAll(Su0Bf[patchi], facei)
198  {
199  Su0Bf[patchi][facei] =
200  Su0pTphi
201  (
202  p.boundaryField()[patchi][facei],
203  Tu.boundaryField()[patchi][facei],
204  phi.boundaryField()[patchi][facei],
205  0.0
206  );
207  }
208  }
209 
210  return tSu0;
211 }
212 
213 
216 {
217  if (psiuReactionThermo_.composition().contains("ft"))
218  {
219  const volScalarField& ft = psiuReactionThermo_.composition().Y("ft");
220 
221  return Su0pTphi
222  (
223  psiuReactionThermo_.p(),
224  psiuReactionThermo_.Tu(),
226  (
227  "stoichiometricAirFuelMassRatio", dimless, psiuReactionThermo_
228  )*ft/max(1 - ft, SMALL)
229  );
230  }
231  else
232  {
233  return Su0pTphi
234  (
235  psiuReactionThermo_.p(),
236  psiuReactionThermo_.Tu(),
237  equivalenceRatio_
238  );
239  }
240 }
241 
242 
243 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::laminarFlameSpeedModels::defineTypeNameAndDebug
defineTypeNameAndDebug(constant, 0)
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::psiuReactionThermo
Foam::psiuReactionThermo.
Definition: psiuReactionThermo.H:55
Foam::laminarFlameSpeedModels::Gulders::~Gulders
virtual ~Gulders()
Destructor.
Definition: Gulders.C:71
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::laminarFlameSpeedModels::Gulders::operator()
tmp< volScalarField > operator()() const
Return the laminar flame speed [m/s].
Definition: Gulders.C:215
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::laminarFlameSpeedModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(laminarFlameSpeed, constant, dictionary)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:749
Foam::laminarFlameSpeed
Abstract class for laminar flame speed.
Definition: laminarFlameSpeed.H:60
Gulders.H
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123