LangmuirHinshelwoodReactionRateI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2012 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  Istream& is
56 )
57 :
58  co_(st["CO"]),
59  c3h6_(st["C3H6"]),
60  no_(st["NO"])
61 {
62  is.readBegin("LangmuirHinshelwoodReactionRate(Istream&)");
63 
64  for (int i=0; i<n_; i++)
65  {
66  is >> A_[i] >> Ta_[i];
67  }
68 
69  is.readEnd("LangmuirHinshelwoodReactionRate(Istream&)");
70 }
71 
72 
74 (
75  const speciesTable& st,
76  const dictionary& dict
77 )
78 :
79  co_(st["CO"]),
80  c3h6_(st["C3H6"]),
81  no_(st["NO"])
82 {
83  // read (A, Ta) pairs
84  FixedList<Tuple2<scalar, scalar>, n_> coeffs(dict.lookup("coeffs"));
85 
86  forAll(coeffs, i)
87  {
88  A_[i] = coeffs[i].first();
89  Ta_[i] = coeffs[i].second();
90  }
91 }
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
96 inline Foam::scalar Foam::LangmuirHinshelwoodReactionRate::operator()
97 (
98  const scalar p,
99  const scalar T,
100  const scalarField& c
101 ) const
102 {
103  return A_[0]*exp(-Ta_[0]/T)/
104  (
105  T
106  *sqr(1 + A_[1]*exp(-Ta_[1]/T)*c[co_] + A_[2]*exp(-Ta_[2]/T)*c[c3h6_])
107  *(1 + A_[3]*exp(-Ta_[3]/T)*sqr(c[co_])*sqr(c[c3h6_]))
108  *(1 + A_[4]*exp(-Ta_[4]/T)*pow(c[no_], 0.7))
109  );
110 }
111 
112 
114 {
116 
117  forAll(coeffs, i)
118  {
119  coeffs[i].first() = A_[i];
120  coeffs[i].second() = Ta_[i];
121  }
122 
123  os.writeKeyword("coeffs") << coeffs << nl;
124 }
125 
126 
127 inline Foam::Ostream& Foam::operator<<
128 (
129  Ostream& os,
131 )
132 {
133  os << token::BEGIN_LIST;
134 
135  for (int i=0; i<LangmuirHinshelwoodReactionRate::n_; i++)
136  {
137  os << token::SPACE
138  << '(' << lhrr.A_[i] << token::SPACE << lhrr.Ta_[i] << ')';
139  }
140 
141  os << token::END_LIST;
142 
143  return os;
144 }
145 
146 
147 // ************************************************************************* //
Istream & readBegin(const char *funcName)
Definition: Istream.C:86
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:53
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensionedSymmTensor sqr(const dimensionedVector &dv)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void write(Ostream &os) const
Write to stream.
Istream & readEnd(const char *funcName)
Definition: Istream.C:103
dimensionedScalar exp(const dimensionedScalar &ds)
T & first()
Return the first element of the list.
Definition: FixedListI.H:207
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
const volScalarField & T
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.
volScalarField & p
LangmuirHinshelwoodReactionRate(const scalar A[], const scalar Ta[], const label co, const label c3h6, const label no)
Construct from components.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451