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-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 // * * * * * * * * * * * * * * * * 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_(readScalar(dict.lookup("A"))),
50  beta_(readScalar(dict.lookup("beta"))),
51  Ta_(readScalar(dict.lookup("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
64 {
65  scalar lta = A_;
66 
67  if (mag(beta_) > vSmall)
68  {
69  lta *= pow(T, beta_);
70  }
71 
72  scalar expArg = 0;
73 
74  if (mag(Ta_) > vSmall)
75  {
76  expArg -= Ta_/T;
77  }
78 
79  scalar lnT = log(T);
80 
81  for (int n=0; n<nb_; n++)
82  {
83  expArg += b_[n]*pow(lnT, n);
84  }
85 
86  lta *= exp(expArg);
87 
88  return lta;
89 }
90 
91 
92 inline Foam::scalar Foam::JanevReactionRate::ddT
93 (
94  const scalar p,
95  const scalar T,
96  const scalarField&
97 ) const
98 {
99  scalar lta = A_;
100 
101  if (mag(beta_) > vSmall)
102  {
103  lta *= pow(T, beta_);
104  }
105 
106  scalar expArg = 0;
107 
108  if (mag(Ta_) > vSmall)
109  {
110  expArg -= Ta_/T;
111  }
112 
113  scalar lnT = log(T);
114 
115  for (int n=0; n<nb_; n++)
116  {
117  expArg += b_[n]*pow(lnT, n);
118  }
119 
120  scalar deriv = b_[1];
121 
122  for (int n=2; n<nb_; n++)
123  {
124  deriv += n*b_[n]*pow(lnT, n-1);
125  }
126 
127  lta *= exp(expArg);
128 
129  return lta*(beta_+Ta_/T+deriv)/T;
130 }
131 
132 
135 {
136  return NullObjectRef<List<Tuple2<label, scalar>>>();
137 }
138 
139 
141 (
142  const scalar p,
143  const scalar T,
144  const scalarField& c,
145  scalarField& dcidc
146 ) const
147 {}
148 
149 
150 inline Foam::scalar Foam::JanevReactionRate::dcidT
151 (
152  const scalar p,
153  const scalar T,
154  const scalarField& c
155 ) const
156 {
157  return 0;
158 }
159 
160 
162 {
163  os.writeKeyword("A") << A_ << nl;
164  os.writeKeyword("beta") << beta_ << nl;
165  os.writeKeyword("Ta") << Ta_ << nl;
166  os.writeKeyword("b") << b_ << nl;
167 }
168 
169 
170 inline Foam::Ostream& Foam::operator<<
171 (
172  Ostream& os,
173  const JanevReactionRate& jrr
174 )
175 {
176  jrr.write(os);
177  return os;
178 }
179 
180 
181 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
scalar ddT(const scalar p, const scalar T, const scalarField &c) const
scalar dcidT(const scalar p, const scalar T, const scalarField &c) const
Temperature derivative of the pressure dependent term.
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)
void dcidc(const scalar p, const scalar T, const scalarField &c, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
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:265
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.
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
virtual Ostream & write(const token &)=0
Write next token to stream.
volScalarField & p
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:583