fluxLimitedLangmuirHinshelwoodReactionRateI.H
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) 2019-2021 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 "volFields.H"
27 #include "mathematicalConstants.H"
28 #include "thermodynamicConstants.H"
29 
30 using namespace Foam::constant::mathematical;
31 using namespace Foam::constant::thermodynamic;
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 inline Foam::scalar Foam::fluxLimitedLangmuirHinshelwoodReactionRate::Av
36 (
37  const label li
38 ) const
39 {
40  if (AvUniform_)
41  {
42  return Av_;
43  }
44  else
45  {
46  return AvField_[li];
47  }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const speciesTable& st,
57  const objectRegistry& ob,
58  const dictionary& dict
59 )
60 :
61  additionalAdsorbableSpecieNames_
62  (
63  dict.lookupOrDefault
64  (
65  "additionalAdsorbableSpecies",
66  wordList()
67  )
68  ),
69  a_(dict.lookupOrDefault("a", scalar(1))),
70  A_(dict.lookup("A")),
71  Ta_(dict.lookup("Ta")),
72  limited_(dict.found("s")),
73  AvUniform_
74  (
75  limited_
76  ? dict.lookup("Av")[0].isNumber()
77  : true
78  ),
79  Av_
80  (
81  limited_ && AvUniform_
82  ? dict.lookup<scalar>("Av")
83  : 0
84  ),
85  AvName_
86  (
87  limited_ && !AvUniform_
88  ? dict.lookup("Av")
89  : word::null
90  ),
91  AvField_
92  (
93  limited_ && !AvUniform_
96  )
97 {
101  (
102  IStringStream(dict.lookup("reaction"))(),
103  st,
104  lhs,
105  rhs
106  );
107 
108  nReactants_ = lhs.size();
109 
110  ra_.setSize(nReactants_ + additionalAdsorbableSpecieNames_.size());
111  nu_.setSize(nReactants_);
112  exp_.setSize(nReactants_);
113 
114  forAll(lhs, i)
115  {
116  ra_[i] = lhs[i].index;
117  nu_[i] = lhs[i].stoichCoeff;
118  exp_[i] = lhs[i].exponent;
119  }
120 
121  forAll(additionalAdsorbableSpecieNames_, i)
122  {
123  ra_[nReactants_ + i] = st[additionalAdsorbableSpecieNames_[i]];
124  }
125 
126  beta_ =
127  dict.lookupOrDefault
128  (
129  "beta",
130  scalarList
131  (
132  1 + nReactants_ + additionalAdsorbableSpecieNames_.size(),
133  scalar(0)
134  )
135  );
136 
137  m_ =
138  dict.lookupOrDefault
139  (
140  "m",
141  scalarList
142  (
143  1 + nReactants_ + additionalAdsorbableSpecieNames_.size(),
144  scalar(1)
145  )
146  );
147 
148  if (!dict.found("m"))
149  {
150  m_[0] = 2;
151  }
152 
153  s_ = dict.lookupOrDefault("s", scalarList(nReactants_, scalar(0)));
154 
155  W_ =
156  limited_
157  ? scalarList(dict.lookup("W"))
158  : scalarList(nReactants_, scalar(0));
159 
160  const scalar nCoeffs =
161  1 + nReactants_ + additionalAdsorbableSpecieNames_.size();
162 
163  if (A_.size() != nCoeffs)
164  {
166  << "Number of A coefficients != " << nCoeffs
167  << exit(FatalIOError);
168  }
169 
170  if (Ta_.size() != nCoeffs)
171  {
173  << "Number of Ta coefficients != " << nCoeffs
174  << exit(FatalIOError);
175  }
176 
177  if (beta_.size() != nCoeffs)
178  {
180  << "Number of beta coefficients != " << nCoeffs
181  << exit(FatalIOError);
182  }
183 
184  if (m_.size() != nCoeffs)
185  {
187  << "Number of m coefficients != " << nCoeffs
188  << exit(FatalIOError);
189  }
190 
191  if (s_.size() != nReactants_)
192  {
194  << "Number of s coefficients != "
195  << nReactants_ << exit(FatalIOError);
196  }
197 
198  if (W_.size() != nReactants_)
199  {
201  << "Number of W coefficients != "
202  << nReactants_ << exit(FatalIOError);
203  }
204 }
205 
206 
207 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
208 
209 inline void
211 {}
212 
213 
214 inline void
216 {}
217 
218 
219 inline Foam::scalar
220 Foam::fluxLimitedLangmuirHinshelwoodReactionRate::operator()
221 (
222  const scalar p,
223  const scalar T,
224  const scalarField& c,
225  const label li
226 ) const
227 {
228  scalar sumkc = 0;
229  forAll(ra_, i)
230  {
231  const label ip1 = i + 1;
232 
233  const scalar kc =
234  A_[ip1]*pow(T, beta_[ip1])*exp(-Ta_[ip1]/T)
235  *pow(c[ra_[i]], m_[ip1]);
236 
237  sumkc += kc;
238  }
239 
240  const scalar TaByT0 = Ta_[0]/T;
241  const scalar k0 = A_[0]*pow(T, beta_[0])*exp(-TaByT0);
242 
243  scalar r = k0/pow(a_ + sumkc, m_[0]);
244 
245  if (limited_)
246  {
247  scalar rc = 1;
248  forAll(ra_, i)
249  {
250  rc *= pow(c[ra_[i]], exp_[i]);
251  }
252 
253  const scalar Av = this->Av(li);
254 
255  if (rc > rootSmall)
256  {
257  forAll(s_, i)
258  {
259  r = min
260  (
261  r,
262  (Av*s_[i]*c[ra_[i]]/(nu_[i]*rc))*sqrt(RR*T/(twoPi*W_[i]))
263  );
264  }
265  }
266  }
267 
268  return r;
269 }
270 
271 
273 (
274  const scalar p,
275  const scalar T,
276  const scalarField& c,
277  const label li
278 ) const
279 {
280  scalar sumkc = 0;
281  scalar sumBetaKc = 0;
282  forAll(ra_, i)
283  {
284  const label ip1 = i + 1;
285 
286  const scalar kc =
287  A_[ip1]*pow(T, beta_[ip1])*exp(-Ta_[ip1]/T)
288  *pow(c[ra_[i]], m_[ip1]);
289 
290  sumkc += kc;
291  sumBetaKc += (beta_[ip1] + Ta_[ip1]/T)*kc;
292  }
293 
294  const scalar TaByT0 = Ta_[0]/T;
295  const scalar k0 = A_[0]*pow(T, beta_[0])*exp(-TaByT0);
296 
297  scalar ddT =
298  ((beta_[0] + TaByT0)*k0 - m_[0]*k0*sumBetaKc/(a_ + sumkc))
299  /(pow(a_ + sumkc, m_[0])*T);
300 
301  if (limited_)
302  {
303  scalar rc = 1;
304  forAll(ra_, i)
305  {
306  rc *= pow(c[ra_[i]], exp_[i]);
307  }
308 
309  scalar r = k0/pow(a_ + sumkc, m_[0]);
310 
311  if (rc > rootSmall)
312  {
313  label l = -1;
314 
315  const scalar Av = this->Av(li);
316 
317  forAll(s_, i)
318  {
319  const scalar rl =
320  (Av*s_[i]*c[ra_[i]]/(nu_[i]*rc))*sqrt(RR*T/(twoPi*W_[i]));
321 
322  if (rl < r)
323  {
324  l = i;
325  r = rl;
326  }
327  }
328 
329  if (l != -1)
330  {
331  ddT =
332  (Av*s_[l]*c[ra_[l]]/(nu_[l]*rc))
333  *0.5*sqrt(RR/(twoPi*W_[l]*T));
334  }
335  }
336  }
337 
338  return ddT;
339 }
340 
341 
344 {
345  return NullObjectRef<List<Tuple2<label, scalar>>>();
346 }
347 
348 
350 (
351  const scalar p,
352  const scalar T,
353  const scalarField& c,
354  const label li,
355  scalarField& dcidc
356 ) const
357 {}
358 
359 
361 (
362  const scalar p,
363  const scalar T,
364  const scalarField& c,
365  const label li
366 ) const
367 {
368  return 0;
369 }
370 
371 
373 (
374  Ostream& os
375 ) const
376 {
377  writeEntry
378  (
379  os,
380  "additionalAdsorbableSpecies",
381  additionalAdsorbableSpecieNames_
382  );
383  writeEntry(os, "a", a_);
384  writeEntry(os, "A", A_);
385  writeEntry(os, "Ta", Ta_);
386  writeEntry(os, "beta", beta_);
387  writeEntry(os, "m", m_);
388  if (limited_)
389  {
390  writeEntry(os, "s", s_);
391  writeEntry(os, "W", W_);
392  }
393 }
394 
395 
396 inline Foam::Ostream& Foam::operator<<
397 (
398  Ostream& os,
400 )
401 {
402  lhrr.write(os);
403  return os;
404 }
405 
406 
407 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
dimensionedScalar sqrt(const dimensionedScalar &ds)
Thermodynamic scalar constants.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const dimensionedScalar c
Speed of light in a vacuum.
fluxLimitedLangmuirHinshelwoodReactionRate(const speciesTable &species, const objectRegistry &ob, const dictionary &dict)
Construct from dictionary.
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
dimensionedScalar exp(const dimensionedScalar &ds)
mathematical constants.
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
static const word null
An empty word.
Definition: word.H:77
const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
Langmuir-Hinshelwood reaction rate for gaseous reactions on surfaces including the optional flux limi...
void dcidc(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
const scalar twoPi(2 *pi)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
scalar dcidT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of the pressure dependent term.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
static void setLRhs(Istream &, const speciesTable &, List< specieCoeffs > &lhs, List< specieCoeffs > &rhs)
Construct the left- and right-hand-side reaction coefficients.
Definition: specieCoeffs.C:95
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:335
A wordList with hashed indices for faster lookup by name.
Input from memory buffer stream.
Definition: IStringStream.H:49
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
volScalarField & p
Registry of regIOobjects.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844
IOerror FatalIOError