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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "FixedList.H"
29#include "Tuple2.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
34(
35 const scalar A[],
36 const scalar Ta[],
37 const label co,
38 const label c3h6,
39 const label no
40)
41:
42 co_(co),
43 c3h6_(c3h6),
44 no_(no)
45{
46 for (int i=0; i<n_; i++)
47 {
48 A_[i] = A[i];
49 Ta_[i] = Ta[i];
50 }
51}
52
53
55(
56 const speciesTable& species,
57 const dictionary& dict
58)
59:
60 co_(species.find("CO")),
61 c3h6_(species.find("C3H6")),
62 no_(species.find("NO"))
63{
64 // read (A, Ta) pairs
65 FixedList<Tuple2<scalar, scalar>, n_> coeffs(dict.lookup("coeffs"));
66
67 forAll(coeffs, i)
68 {
69 A_[i] = coeffs[i].first();
70 Ta_[i] = coeffs[i].second();
71 }
72}
73
74
75// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76
78(
79 const scalar p,
80 const scalar T,
81 const scalarField& c
82) const
83{
84 return A_[0]*exp(-Ta_[0]/T)/
85 (
86 T
87 *sqr(1 + A_[1]*exp(-Ta_[1]/T)*c[co_] + A_[2]*exp(-Ta_[2]/T)*c[c3h6_])
88 *(1 + A_[3]*exp(-Ta_[3]/T)*sqr(c[co_])*sqr(c[c3h6_]))
89 *(1 + A_[4]*exp(-Ta_[4]/T)*pow(c[no_], 0.7))
90 );
91}
92
93
95{
97
98 forAll(coeffs, i)
99 {
100 coeffs[i].first() = A_[i];
101 coeffs[i].second() = Ta_[i];
102 }
103
104 os.writeKeyword("coeffs") << coeffs << nl;
105}
106
107
108inline Foam::Ostream& Foam::operator<<
109(
110 Ostream& os,
112)
113{
114 lhrr.write(os);
115 return os;
116}
117
118
119// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: FixedList.H:81
T & first() noexcept
The first element of the list, position [0].
Definition: FixedListI.H:207
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
virtual Ostream & writeKeyword(const keyType &kw)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:57
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:386
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
virtual bool write()
Write the output fields.
A wordList with hashed named lookup, which can be faster in some situations than using the normal lis...
volScalarField & p
const volScalarField & T
OBJstream os(runTime.globalPath()/outputName)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333