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-2018 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 scalar A[],
34  const scalar Ta[],
35  const label co,
36  const label c3h6,
37  const label no
38 )
39 :
40  co_(co),
41  c3h6_(c3h6),
42  no_(no)
43 {
44  for (int i=0; i<n_; i++)
45  {
46  A_[i] = A[i];
47  Ta_[i] = Ta[i];
48  }
49 }
50 
51 
53 (
54  const speciesTable& st,
55  const dictionary& dict
56 )
57 :
58  co_(st["CO"]),
59  c3h6_(st["C3H6"]),
60  no_(st["NO"])
61 {
62  // read (A, Ta) pairs
63  FixedList<Tuple2<scalar, scalar>, n_> coeffs(dict.lookup("coeffs"));
64 
65  forAll(coeffs, i)
66  {
67  A_[i] = coeffs[i].first();
68  Ta_[i] = coeffs[i].second();
69  }
70 }
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 inline Foam::scalar Foam::LangmuirHinshelwoodReactionRate::operator()
76 (
77  const scalar p,
78  const scalar T,
79  const scalarField& c
80 ) const
81 {
82  return A_[0]*exp(-Ta_[0]/T)/
83  (
84  T
85  *sqr(1 + A_[1]*exp(-Ta_[1]/T)*c[co_] + A_[2]*exp(-Ta_[2]/T)*c[c3h6_])
86  *(1 + A_[3]*exp(-Ta_[3]/T)*sqr(c[co_])*sqr(c[c3h6_]))
87  *(1 + A_[4]*exp(-Ta_[4]/T)*pow(c[no_], 0.7))
88  );
89 }
90 
91 
93 (
94  const scalar p,
95  const scalar T,
96  const scalarField& c
97 ) const
98 {
99  scalar den =
100  (
101  T
102  *sqr(1 + A_[1]*exp(-Ta_[1]/T)*c[co_] + A_[2]*exp(-Ta_[2]/T)*c[c3h6_])
103  *(1 + A_[3]*exp(-Ta_[3]/T)*sqr(c[co_])*sqr(c[c3h6_]))
104  *(1 + A_[4]*exp(-Ta_[4]/T)*pow(c[no_], 0.7))
105  );
106  scalar rate = A_[0]*exp(-Ta_[0]/T)/ den;
107 
108  scalar derivDen =
109  (
110  sqr(1 + A_[1]*exp(-Ta_[1]/T)*c[co_] + A_[2]*exp(-Ta_[2]/T)*c[c3h6_])
111  *(1 + A_[3]*exp(-Ta_[3]/T)*sqr(c[co_])*sqr(c[c3h6_]))
112  *(1 + A_[4]*exp(-Ta_[4]/T)*pow(c[no_], 0.7))
113  + 2*T*(1 + A_[1]*exp(-Ta_[1]/T)*c[co_] + A_[2]*exp(-Ta_[2]/T)*c[c3h6_])
114  *
115  (
116  A_[1]*exp(-Ta_[1]/T)*c[co_]*Ta_[1]/sqr(T)
117  + A_[2]*exp(-Ta_[2]/T)*c[c3h6_]*Ta_[2]/sqr(T)
118  )
119  *(1 + A_[3]*exp(-Ta_[3]/T)*sqr(c[co_])*sqr(c[c3h6_]))
120  *(1 + A_[4]*exp(-Ta_[4]/T)*pow(c[no_], 0.7))
121  + T
122  *sqr(1 + A_[1]*exp(-Ta_[1]/T)*c[co_] + A_[2]*exp(-Ta_[2]/T)*c[c3h6_])
123  *(A_[3]*exp(-Ta_[3]/T)*Ta_[3]*sqr(c[co_])*sqr(c[c3h6_])/sqr(T))
124  *(1 + A_[4]*exp(-Ta_[4]/T)*pow(c[no_], 0.7))
125  + T
126  *sqr(1 + A_[1]*exp(-Ta_[1]/T)*c[co_] + A_[2]*exp(-Ta_[2]/T)*c[c3h6_])
127  *(1 + A_[3]*exp(-Ta_[3]/T)*sqr(c[co_])*sqr(c[c3h6_]))
128  *(A_[4]*exp(-Ta_[4]/T)*Ta_[4]*pow(c[no_], 0.7))/sqr(T)
129  );
130 
131  return rate*(Ta_[0]/sqr(T) - derivDen/den);
132 }
133 
134 
137 {
138  return NullObjectRef<List<Tuple2<label, scalar>>>();
139 }
140 
141 
143 (
144  const scalar p,
145  const scalar T,
146  const scalarField& c,
147  scalarField& dcidc
148 ) const
149 {}
150 
151 
153 (
154  const scalar p,
155  const scalar T,
156  const scalarField& c
157 ) const
158 {
159  return 0;
160 }
161 
162 
164 {
166 
167  forAll(coeffs, i)
168  {
169  coeffs[i].first() = A_[i];
170  coeffs[i].second() = Ta_[i];
171  }
172 
173  os.writeKeyword("coeffs") << coeffs << nl;
174 }
175 
176 
177 inline Foam::Ostream& Foam::operator<<
178 (
179  Ostream& os,
181 )
182 {
183  lhrr.write(os);
184  return os;
185 }
186 
187 
188 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
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:137
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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)
T & first()
Return the first element of the list.
Definition: FixedListI.H:227
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)
static const char nl
Definition: Ostream.H:265
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
const volScalarField & T
void write(Ostream &os) const
Write to stream.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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:576
LangmuirHinshelwoodReactionRate(const scalar A[], const scalar Ta[], const label co, const label c3h6, const label no)
Construct from components.