JanevReactionRateI.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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
29 (
30  const scalar A,
31  const scalar beta,
32  const scalar Ta,
34 )
35 :
36  A_(A),
37  beta_(beta),
38  Ta_(Ta),
39  b_(b)
40 {}
41 
42 
44 (
45  const speciesTable&,
46  const dictionary& dict
47 )
48 :
49  A_(dict.lookup<scalar>("A")),
50  beta_(dict.lookup<scalar>("beta")),
51  Ta_(dict.lookup<scalar>("Ta")),
52  b_(dict.lookup("b"))
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
57 
58 inline Foam::scalar Foam::JanevReactionRate::operator()
59 (
60  const scalar p,
61  const scalar T,
62  const scalarField&,
63  const label
64 ) const
65 {
66  scalar lta = A_;
67 
68  if (mag(beta_) > vSmall)
69  {
70  lta *= pow(T, beta_);
71  }
72 
73  scalar expArg = 0;
74 
75  if (mag(Ta_) > vSmall)
76  {
77  expArg -= Ta_/T;
78  }
79 
80  scalar lnT = log(T);
81 
82  for (int n=0; n<nb_; n++)
83  {
84  expArg += b_[n]*pow(lnT, n);
85  }
86 
87  lta *= exp(expArg);
88 
89  return lta;
90 }
91 
92 
93 inline Foam::scalar Foam::JanevReactionRate::ddT
94 (
95  const scalar p,
96  const scalar T,
97  const scalarField&,
98  const label
99 ) const
100 {
101  scalar lta = A_;
102 
103  if (mag(beta_) > vSmall)
104  {
105  lta *= pow(T, beta_);
106  }
107 
108  scalar expArg = 0;
109 
110  if (mag(Ta_) > vSmall)
111  {
112  expArg -= Ta_/T;
113  }
114 
115  scalar lnT = log(T);
116 
117  for (int n=0; n<nb_; n++)
118  {
119  expArg += b_[n]*pow(lnT, n);
120  }
121 
122  scalar deriv = b_[1];
123 
124  for (int n=2; n<nb_; n++)
125  {
126  deriv += n*b_[n]*pow(lnT, n-1);
127  }
128 
129  lta *= exp(expArg);
130 
131  return lta*(beta_+Ta_/T+deriv)/T;
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  const label li,
148  scalarField& dcidc
149 ) const
150 {}
151 
152 
153 inline Foam::scalar Foam::JanevReactionRate::dcidT
154 (
155  const scalar p,
156  const scalar T,
157  const scalarField& c,
158  const label li
159 ) const
160 {
161  return 0;
162 }
163 
164 
166 {
167  writeKeyword(os, "A") << A_ << nl;
168  writeKeyword(os, "beta") << beta_ << nl;
169  writeKeyword(os, "Ta") << Ta_ << nl;
170  writeKeyword(os, "b") << b_ << nl;
171 }
172 
173 
174 inline Foam::Ostream& Foam::operator<<
175 (
176  Ostream& os,
177  const JanevReactionRate& jrr
178 )
179 {
180  jrr.write(os);
181  return os;
182 }
183 
184 
185 // ************************************************************************* //
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
Ostream & writeKeyword(Foam::Ostream &os, const keyType &kw)
Write the keyword to the Ostream with the current level of indentation.
Definition: keyType.C:155
dimensionedScalar log(const dimensionedScalar &ds)
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
const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
Janev, Langer, Evans and Post reaction rate.
dimensionedScalar exp(const dimensionedScalar &ds)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
const volScalarField & T
scalar dcidT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of the pressure dependent term.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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.
A wordList with hashed indices for faster lookup by name.
JanevReactionRate(const scalar A, const scalar beta, const scalar Ta, const FixedList< scalar, nb_ > b)
Construct from components.
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
volScalarField & p
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
void write(Ostream &os) const
Write to stream.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812