35 inline Foam::scalar Foam::fluxLimitedLangmuirHinshelwoodReactionRate::Av
63 additionalAdsorbableSpecieNames_
67 "additionalAdsorbableSpecies",
72 a_(
dict.lookupOrDefault(
"a", scalar(1))),
74 Ta_(
dict.lookup(
"Ta")),
78 AvUniform_(limited_ ?
dict.lookup(
"Av")[0].isNumber() : true),
79 Av_(limited_ && AvUniform_ ?
dict.lookup<scalar>(
"Av") : 0),
80 AvName_(limited_ && !AvUniform_ ?
dict.lookup(
"Av") :
word::null),
98 nReactants_ = lhs.
size();
100 ra_.
setSize(nReactants_ + additionalAdsorbableSpecieNames_.
size());
106 ra_[i] = lhs[i].index;
107 nu_[i] = lhs[i].stoichCoeff;
108 exp_[i] = lhs[i].exponent;
111 forAll(additionalAdsorbableSpecieNames_, i)
113 ra_[nReactants_ + i] = st[additionalAdsorbableSpecieNames_[i]];
122 1 + nReactants_ + additionalAdsorbableSpecieNames_.
size(),
133 1 + nReactants_ + additionalAdsorbableSpecieNames_.
size(),
138 if (!
dict.found(
"m"))
143 s_ =
dict.lookupOrDefault(
"s",
scalarList(nReactants_, scalar(0)));
150 const scalar nCoeffs =
151 1 + nReactants_ + additionalAdsorbableSpecieNames_.
size();
153 if (A_.
size() != nCoeffs)
156 <<
"Number of A coefficients != " << nCoeffs
160 if (Ta_.
size() != nCoeffs)
163 <<
"Number of Ta coefficients != " << nCoeffs
167 if (beta_.
size() != nCoeffs)
170 <<
"Number of beta coefficients != " << nCoeffs
174 if (m_.
size() != nCoeffs)
177 <<
"Number of m coefficients != " << nCoeffs
181 if (s_.
size() != nReactants_)
184 <<
"Number of s coefficients != "
188 if (W_.
size() != nReactants_)
191 <<
"Number of W coefficients != "
220 Foam::fluxLimitedLangmuirHinshelwoodReactionRate::operator()
231 const label ip1 = i + 1;
234 A_[ip1]*
pow(
T, beta_[ip1])*
exp(-Ta_[ip1]/
T)
235 *
pow(
c[ra_[i]], m_[ip1]);
240 const scalar TaByT0 = Ta_[0]/
T;
241 const scalar k0 = A_[0]*
pow(
T, beta_[0])*
exp(-TaByT0);
243 scalar r = k0/
max(
pow(a_ + sumkc, m_[0]), rootVSmall);
250 rc *=
pow(
c[ra_[i]], exp_[i]);
266 return this->Av(li)*r;
279 scalar sumBetaKc = 0;
282 const label ip1 = i + 1;
285 A_[ip1]*
pow(
T, beta_[ip1])*
exp(-Ta_[ip1]/
T)
286 *
pow(
c[ra_[i]], m_[ip1]);
289 sumBetaKc += (beta_[ip1] + Ta_[ip1]/
T)*kc;
292 const scalar TaByT0 = Ta_[0]/
T;
293 const scalar k0 = A_[0]*
pow(
T, beta_[0])*
exp(-TaByT0);
296 ((beta_[0] + TaByT0)*k0*(a_ + sumkc) - m_[0]*k0*sumBetaKc)
297 /(
max(
pow(a_ + sumkc, m_[0] + 1), rootVSmall)*
T);
304 rc *=
pow(
c[ra_[i]], exp_[i]);
307 scalar r = k0/
pow(a_ + sumkc, m_[0]);
328 (s_[l]*
c[ra_[l]]/(nu_[l]*rc))
334 return this->Av(li)*ddT;
365 "additionalAdsorbableSpecies",
366 additionalAdsorbableSpecieNames_
#define forAll(list, i)
Loop across all elements in list.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Input from memory buffer stream.
void size(const label)
Override size to be inconsistent with allocated storage.
void setSize(const label)
Reset size of List.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
virtual Ostream & write(const char)=0
Write character.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Dimension set for the base types.
Langmuir-Hinshelwood reaction rate for gaseous reactions on surfaces including the optional flux limi...
void preEvaluate() const
Pre-evaluation hook.
void postEvaluate() const
Post-evaluation hook.
void write(Ostream &os) const
Write to stream.
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
fluxLimitedLangmuirHinshelwoodReactionRate(const speciesTable &species, const objectRegistry &ob, const dimensionSet &dims, const dictionary &dict)
Construct from dictionary.
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.
A class for managing temporary objects.
A class for handling words, derived from string.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
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)
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.
List< scalar > scalarList
A List of scalars.
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)
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)