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-2022 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/max(pow(a_ + sumkc, m_[0]), rootVSmall);
244 
245  if (limited_)
246  {
247  scalar rc = 1;
248  forAll(exp_, i)
249  {
250  rc *= pow(c[ra_[i]], exp_[i]);
251  }
252 
253  if (rc > rootSmall)
254  {
255  forAll(s_, i)
256  {
257  r = min
258  (
259  r,
260  (s_[i]*c[ra_[i]]/(nu_[i]*rc))*sqrt(RR*T/(twoPi*W_[i]))
261  );
262  }
263  }
264  }
265 
266  return this->Av(li)*r;
267 }
268 
269 
271 (
272  const scalar p,
273  const scalar T,
274  const scalarField& c,
275  const label li
276 ) const
277 {
278  scalar sumkc = 0;
279  scalar sumBetaKc = 0;
280  forAll(ra_, i)
281  {
282  const label ip1 = i + 1;
283 
284  const scalar kc =
285  A_[ip1]*pow(T, beta_[ip1])*exp(-Ta_[ip1]/T)
286  *pow(c[ra_[i]], m_[ip1]);
287 
288  sumkc += kc;
289  sumBetaKc += (beta_[ip1] + Ta_[ip1]/T)*kc;
290  }
291 
292  const scalar TaByT0 = Ta_[0]/T;
293  const scalar k0 = A_[0]*pow(T, beta_[0])*exp(-TaByT0);
294 
295  scalar ddT =
296  ((beta_[0] + TaByT0)*k0*(a_ + sumkc) - m_[0]*k0*sumBetaKc)
297  /(max(pow(a_ + sumkc, m_[0] + 1), rootVSmall)*T);
298 
299  if (limited_)
300  {
301  scalar rc = 1;
302  forAll(exp_, i)
303  {
304  rc *= pow(c[ra_[i]], exp_[i]);
305  }
306 
307  scalar r = k0/pow(a_ + sumkc, m_[0]);
308 
309  if (rc > rootSmall)
310  {
311  label l = -1;
312 
313  forAll(s_, i)
314  {
315  const scalar rl =
316  (s_[i]*c[ra_[i]]/(nu_[i]*rc))*sqrt(RR*T/(twoPi*W_[i]));
317 
318  if (rl < r)
319  {
320  l = i;
321  r = rl;
322  }
323  }
324 
325  if (l != -1)
326  {
327  ddT =
328  (s_[l]*c[ra_[l]]/(nu_[l]*rc))
329  *0.5*sqrt(RR/(twoPi*W_[l]*T));
330  }
331  }
332  }
333 
334  return this->Av(li)*ddT;
335 }
336 
337 
339 {
340  return false;
341 }
342 
343 
345 (
346  const scalar p,
347  const scalar T,
348  const scalarField& c,
349  const label li,
350  scalarField& ddc
351 ) const
352 {
353  ddc = 0;
354 }
355 
356 
358 (
359  Ostream& os
360 ) const
361 {
362  writeEntry
363  (
364  os,
365  "additionalAdsorbableSpecies",
366  additionalAdsorbableSpecieNames_
367  );
368  writeEntry(os, "a", a_);
369  writeEntry(os, "A", A_);
370  writeEntry(os, "Ta", Ta_);
371  writeEntry(os, "beta", beta_);
372  writeEntry(os, "m", m_);
373  if (limited_)
374  {
375  writeEntry(os, "s", s_);
376  writeEntry(os, "W", W_);
377  }
378 }
379 
380 
381 inline Foam::Ostream& Foam::operator<<
382 (
383  Ostream& os,
385 )
386 {
387  lhrr.write(os);
388  return os;
389 }
390 
391 
392 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual Ostream & write(const char)=0
Write character.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
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
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Langmuir-Hinshelwood reaction rate for gaseous reactions on surfaces including the optional flux limi...
const scalar twoPi(2 *pi)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
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
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:318
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...
void ddc(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &ddc) const
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:864
IOerror FatalIOError