LangmuirHinshelwoodReactionRateI.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) 2011-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 "FixedList.H"
27 #include "Tuple2.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const speciesTable& st,
34  const dictionary& dict
35 )
36 :
37  reactantNames_(dict.lookup("reactants")),
38  a_(dict.lookupOrDefault<scalar>("a", 1)),
39  A_(dict.lookup("A")),
40  Ta_(dict.lookup("Ta")),
41  beta_
42  (
44  (
45  "beta",
46  FixedList<scalar, 3>({0, 0, 0})
47  )
48  ),
49  m_
50  (
52  (
53  "m",
54  FixedList<scalar, 3>({2, 1, 1})
55  )
56  )
57 {
58  forAll(reactantNames_, i)
59  {
60  r_[i] = st[reactantNames_[i]];
61  }
62 }
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
68 {}
69 
70 
72 {}
73 
74 
75 inline Foam::scalar Foam::LangmuirHinshelwoodReactionRate::operator()
76 (
77  const scalar p,
78  const scalar T,
79  const scalarField& c,
80  const label li
81 ) const
82 {
83  const scalar c1m = pow(c[r_[0]], m_[1]);
84  const scalar c2m = pow(c[r_[1]], m_[2]);
85 
86  const scalar TaByT0 = Ta_[0]/T;
87  const scalar TaByT1 = Ta_[1]/T;
88  const scalar TaByT2 = Ta_[2]/T;
89 
90  const scalar k0 = A_[0]*pow(T, beta_[0])*exp(-TaByT0);
91  const scalar k1 = A_[1]*pow(T, beta_[1])*exp(-TaByT1);
92  const scalar k2 = A_[2]*pow(T, beta_[2])*exp(-TaByT2);
93 
94  return k0/pow(a_ + k1*c1m + k2*c2m, m_[0]);
95 }
96 
97 
99 (
100  const scalar p,
101  const scalar T,
102  const scalarField& c,
103  const label li
104 ) const
105 {
106  const scalar c1m = pow(c[r_[0]], m_[1]);
107  const scalar c2m = pow(c[r_[1]], m_[2]);
108 
109  const scalar TaByT0 = Ta_[0]/T;
110  const scalar TaByT1 = Ta_[1]/T;
111  const scalar TaByT2 = Ta_[2]/T;
112 
113  const scalar k0 = A_[0]*pow(T, beta_[0])*exp(-TaByT0);
114  const scalar k1 = A_[1]*pow(T, beta_[1])*exp(-TaByT1);
115  const scalar k2 = A_[2]*pow(T, beta_[2])*exp(-TaByT2);
116 
117  return
118  (
119  (beta_[0] + TaByT0)*k0
120  - m_[0]*k0*((beta_[1] + TaByT1)*k1*c1m + (beta_[2] + TaByT2)*k2*c2m)
121  /(a_ + k1*c1m + k2*c2m)
122  )/(pow(a_ + k1*c1m + k2*c2m, m_[0])*T);
123 }
124 
125 
128 {
129  return NullObjectRef<List<Tuple2<label, scalar>>>();
130 }
131 
132 
134 (
135  const scalar p,
136  const scalar T,
137  const scalarField& c,
138  const label li,
139  scalarField& dcidc
140 ) const
141 {}
142 
143 
145 (
146  const scalar p,
147  const scalar T,
148  const scalarField& c,
149  const label li
150 ) const
151 {
152  return 0;
153 }
154 
155 
157 {
158  writeEntry(os, "reactants", reactantNames_);
159  writeEntry(os, "a", a_);
160  writeEntry(os, "A", A_);
161  writeEntry(os, "Ta", Ta_);
162  writeEntry(os, "beta", beta_);
163  writeEntry(os, "m", m_);
164 }
165 
166 
167 inline Foam::Ostream& Foam::operator<<
168 (
169  Ostream& os,
171 )
172 {
173  lhrr.write(os);
174  return os;
175 }
176 
177 
178 // ************************************************************************* //
#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
Langmuir-Hinshelwood reaction rate for gaseous reactions on surfaces.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar exp(const dimensionedScalar &ds)
void postEvaluate() const
Post-evaluation hook.
LangmuirHinshelwoodReactionRate(const speciesTable &species, const dictionary &dict)
Construct from dictionary.
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
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
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.
void write(Ostream &os) const
Write to stream.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
A wordList with hashed indices for faster lookup by name.
volScalarField & p
scalar dcidT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of the pressure dependent term.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844