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-2023 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 tAv_()[li];
47  }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const speciesTable& st,
57  const objectRegistry& ob,
58  const dictionary& dict
59 )
60 :
61  nReactants_(),
62  additionalAdsorbableSpecieNames_
63  (
64  dict.lookupOrDefault
65  (
66  "additionalAdsorbableSpecies",
67  wordList()
68  )
69  ),
70  ra_(),
71  a_(dict.lookupOrDefault("a", scalar(1))),
72  A_(dict.lookup("A")),
73  Ta_(dict.lookup("Ta")),
74  beta_(),
75  m_(),
76  limited_(dict.found("s")),
77  AvUniform_(limited_ ? dict.lookup("Av")[0].isNumber() : true),
78  Av_(limited_ && AvUniform_ ? dict.lookup<scalar>("Av") : 0),
79  AvName_(limited_ && !AvUniform_ ? dict.lookup("Av") : word::null),
80  nu_(),
81  exp_(),
82  s_(),
83  W_(),
84  ob_(ob),
85  tAv_(nullptr)
86 {
90  (
91  IStringStream(dict.lookup("reaction"))(),
92  st,
93  lhs,
94  rhs
95  );
96 
97  nReactants_ = lhs.size();
98 
99  ra_.setSize(nReactants_ + additionalAdsorbableSpecieNames_.size());
100  nu_.setSize(nReactants_);
101  exp_.setSize(nReactants_);
102 
103  forAll(lhs, i)
104  {
105  ra_[i] = lhs[i].index;
106  nu_[i] = lhs[i].stoichCoeff;
107  exp_[i] = lhs[i].exponent;
108  }
109 
110  forAll(additionalAdsorbableSpecieNames_, i)
111  {
112  ra_[nReactants_ + i] = st[additionalAdsorbableSpecieNames_[i]];
113  }
114 
115  beta_ =
116  dict.lookupOrDefault
117  (
118  "beta",
119  scalarList
120  (
121  1 + nReactants_ + additionalAdsorbableSpecieNames_.size(),
122  scalar(0)
123  )
124  );
125 
126  m_ =
127  dict.lookupOrDefault
128  (
129  "m",
130  scalarList
131  (
132  1 + nReactants_ + additionalAdsorbableSpecieNames_.size(),
133  scalar(1)
134  )
135  );
136 
137  if (!dict.found("m"))
138  {
139  m_[0] = 2;
140  }
141 
142  s_ = dict.lookupOrDefault("s", scalarList(nReactants_, scalar(0)));
143 
144  W_ =
145  limited_
146  ? scalarList(dict.lookup("W"))
147  : scalarList(nReactants_, scalar(0));
148 
149  const scalar nCoeffs =
150  1 + nReactants_ + additionalAdsorbableSpecieNames_.size();
151 
152  if (A_.size() != nCoeffs)
153  {
155  << "Number of A coefficients != " << nCoeffs
156  << exit(FatalIOError);
157  }
158 
159  if (Ta_.size() != nCoeffs)
160  {
162  << "Number of Ta coefficients != " << nCoeffs
163  << exit(FatalIOError);
164  }
165 
166  if (beta_.size() != nCoeffs)
167  {
169  << "Number of beta coefficients != " << nCoeffs
170  << exit(FatalIOError);
171  }
172 
173  if (m_.size() != nCoeffs)
174  {
176  << "Number of m coefficients != " << nCoeffs
177  << exit(FatalIOError);
178  }
179 
180  if (s_.size() != nReactants_)
181  {
183  << "Number of s coefficients != "
184  << nReactants_ << exit(FatalIOError);
185  }
186 
187  if (W_.size() != nReactants_)
188  {
190  << "Number of W coefficients != "
191  << nReactants_ << exit(FatalIOError);
192  }
193 }
194 
195 
196 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
197 
198 inline void
200 {
201  if (!AvUniform_)
202  {
203  tAv_ = ob_.lookupObject<volScalarField::Internal>(AvName_);
204  }
205 }
206 
207 
208 inline void
210 {
211  if (!AvUniform_)
212  {
213  tAv_ = tmp<volScalarField::Internal>(nullptr);
214  }
215 }
216 
217 
218 inline Foam::scalar
219 Foam::fluxLimitedLangmuirHinshelwoodReactionRate::operator()
220 (
221  const scalar p,
222  const scalar T,
223  const scalarField& c,
224  const label li
225 ) const
226 {
227  scalar sumkc = 0;
228  forAll(ra_, i)
229  {
230  const label ip1 = i + 1;
231 
232  const scalar kc =
233  A_[ip1]*pow(T, beta_[ip1])*exp(-Ta_[ip1]/T)
234  *pow(c[ra_[i]], m_[ip1]);
235 
236  sumkc += kc;
237  }
238 
239  const scalar TaByT0 = Ta_[0]/T;
240  const scalar k0 = A_[0]*pow(T, beta_[0])*exp(-TaByT0);
241 
242  scalar r = k0/max(pow(a_ + sumkc, m_[0]), rootVSmall);
243 
244  if (limited_)
245  {
246  scalar rc = 1;
247  forAll(exp_, i)
248  {
249  rc *= pow(c[ra_[i]], exp_[i]);
250  }
251 
252  if (rc > rootSmall)
253  {
254  forAll(s_, i)
255  {
256  r = min
257  (
258  r,
259  (s_[i]*c[ra_[i]]/(nu_[i]*rc))*sqrt(RR*T/(twoPi*W_[i]))
260  );
261  }
262  }
263  }
264 
265  return this->Av(li)*r;
266 }
267 
268 
270 (
271  const scalar p,
272  const scalar T,
273  const scalarField& c,
274  const label li
275 ) const
276 {
277  scalar sumkc = 0;
278  scalar sumBetaKc = 0;
279  forAll(ra_, i)
280  {
281  const label ip1 = i + 1;
282 
283  const scalar kc =
284  A_[ip1]*pow(T, beta_[ip1])*exp(-Ta_[ip1]/T)
285  *pow(c[ra_[i]], m_[ip1]);
286 
287  sumkc += kc;
288  sumBetaKc += (beta_[ip1] + Ta_[ip1]/T)*kc;
289  }
290 
291  const scalar TaByT0 = Ta_[0]/T;
292  const scalar k0 = A_[0]*pow(T, beta_[0])*exp(-TaByT0);
293 
294  scalar ddT =
295  ((beta_[0] + TaByT0)*k0*(a_ + sumkc) - m_[0]*k0*sumBetaKc)
296  /(max(pow(a_ + sumkc, m_[0] + 1), rootVSmall)*T);
297 
298  if (limited_)
299  {
300  scalar rc = 1;
301  forAll(exp_, i)
302  {
303  rc *= pow(c[ra_[i]], exp_[i]);
304  }
305 
306  scalar r = k0/pow(a_ + sumkc, m_[0]);
307 
308  if (rc > rootSmall)
309  {
310  label l = -1;
311 
312  forAll(s_, i)
313  {
314  const scalar rl =
315  (s_[i]*c[ra_[i]]/(nu_[i]*rc))*sqrt(RR*T/(twoPi*W_[i]));
316 
317  if (rl < r)
318  {
319  l = i;
320  r = rl;
321  }
322  }
323 
324  if (l != -1)
325  {
326  ddT =
327  (s_[l]*c[ra_[l]]/(nu_[l]*rc))
328  *0.5*sqrt(RR/(twoPi*W_[l]*T));
329  }
330  }
331  }
332 
333  return this->Av(li)*ddT;
334 }
335 
336 
338 {
339  return false;
340 }
341 
342 
344 (
345  const scalar p,
346  const scalar T,
347  const scalarField& c,
348  const label li,
349  scalarField& ddc
350 ) const
351 {
352  ddc = 0;
353 }
354 
355 
357 (
358  Ostream& os
359 ) const
360 {
361  writeEntry
362  (
363  os,
364  "additionalAdsorbableSpecies",
365  additionalAdsorbableSpecieNames_
366  );
367  writeEntry(os, "a", a_);
368  writeEntry(os, "A", A_);
369  writeEntry(os, "Ta", Ta_);
370  writeEntry(os, "beta", beta_);
371  writeEntry(os, "m", m_);
372  if (limited_)
373  {
374  writeEntry(os, "s", s_);
375  writeEntry(os, "W", W_);
376  }
377 }
378 
379 
380 inline Foam::Ostream& Foam::operator<<
381 (
382  Ostream& os,
384 )
385 {
386  lhrr.write(os);
387  return os;
388 }
389 
390 
391 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Input from memory buffer stream.
Definition: IStringStream.H:52
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual Ostream & write(const char)=0
Write character.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Langmuir-Hinshelwood reaction rate for gaseous reactions on surfaces including the optional flux limi...
fluxLimitedLangmuirHinshelwoodReactionRate(const speciesTable &species, const objectRegistry &ob, const dictionary &dict)
Construct from dictionary.
void ddc(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &ddc) const
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
A wordList with hashed indices for faster lookup by name.
Registry of regIOobjects.
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
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
mathematical constants.
const scalar twoPi(2 *pi)
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
Thermodynamic scalar constants.
const dimensionedScalar c
Speed of light in a vacuum.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar exp(const dimensionedScalar &ds)
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
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p