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-2019 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>({1, 1, 2})
55  )
56  )
57 {
58  forAll(reactantNames_, i)
59  {
60  r_[i] = st[reactantNames_[i]];
61  }
62 }
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
67 inline Foam::scalar Foam::LangmuirHinshelwoodReactionRate::operator()
68 (
69  const scalar p,
70  const scalar T,
71  const scalarField& c
72 ) const
73 {
74  const scalar c0m = pow(c[r_[0]], m_[0]);
75  const scalar c1m = pow(c[r_[1]], m_[1]);
76 
77  const scalar TaByT0 = Ta_[0]/T;
78  const scalar TaByT1 = Ta_[1]/T;
79  const scalar TaByT2 = Ta_[2]/T;
80 
81  const scalar k0 = A_[0]*pow(T, beta_[0])*exp(-TaByT0);
82  const scalar k1 = A_[1]*pow(T, beta_[1])*exp(-TaByT1);
83  const scalar k2 = A_[2]*pow(T, beta_[2])*exp(-TaByT2);
84 
85  return k0/pow(a_ + k1*c0m + k2*c1m, m_[2]);
86 }
87 
88 
90 (
91  const scalar p,
92  const scalar T,
93  const scalarField& c
94 ) const
95 {
96  const scalar c0m = pow(c[r_[0]], m_[0]);
97  const scalar c1m = pow(c[r_[1]], m_[1]);
98 
99  const scalar TaByT0 = Ta_[0]/T;
100  const scalar TaByT1 = Ta_[1]/T;
101  const scalar TaByT2 = Ta_[2]/T;
102 
103  const scalar k0 = A_[0]*pow(T, beta_[0])*exp(-TaByT0);
104  const scalar k1 = A_[1]*pow(T, beta_[1])*exp(-TaByT1);
105  const scalar k2 = A_[2]*pow(T, beta_[2])*exp(-TaByT2);
106 
107  return
108  (
109  (beta_[0] + TaByT0)*k0
110  - m_[2]*k0*((beta_[1] + TaByT1)*k1*c0m + (beta_[2] + TaByT2)*k2*c1m)
111  /(a_ + k1*c0m + k2*c1m)
112  )/(pow(a_ + k1*c0m + k2*c1m, m_[2])*T);
113 }
114 
115 
118 {
119  return NullObjectRef<List<Tuple2<label, scalar>>>();
120 }
121 
122 
124 (
125  const scalar p,
126  const scalar T,
127  const scalarField& c,
128  scalarField& dcidc
129 ) const
130 {}
131 
132 
134 (
135  const scalar p,
136  const scalar T,
137  const scalarField& c
138 ) const
139 {
140  return 0;
141 }
142 
143 
145 {
146  writeEntry(os, "reactants", reactantNames_);
147  writeEntry(os, "a", a_);
148  writeEntry(os, "A", A_);
149  writeEntry(os, "Ta", Ta_);
150  writeEntry(os, "beta", beta_);
151  writeEntry(os, "m", m_);
152 }
153 
154 
155 inline Foam::Ostream& Foam::operator<<
156 (
157  Ostream& os,
159 )
160 {
161  lhrr.write(os);
162  return os;
163 }
164 
165 
166 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Langmuir-Hinshelwood reaction rate for gaseous reactions on surfaces.
scalar ddT(const scalar p, const scalar T, const scalarField &c) const
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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
void dcidc(const scalar p, const scalar T, const scalarField &c, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
dimensionedScalar exp(const dimensionedScalar &ds)
LangmuirHinshelwoodReactionRate(const speciesTable &species, const dictionary &dict)
Construct from dictionary.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
scalar dcidT(const scalar p, const scalar T, const scalarField &c) const
Temperature derivative of the pressure dependent term.
const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const volScalarField & T
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.
const dimensionedScalar c
Speed of light in a vacuum.
virtual Ostream & write(const token &)=0
Write next token to stream.
volScalarField & p
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583