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-2024 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 "JanevReactionRate.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const scalar A,
34  const scalar beta,
35  const scalar Ta,
37 )
38 :
39  beta_(beta),
40  A_(A),
41  Ta_(Ta),
42  b_(b)
43 {}
44 
45 
47 (
48  const speciesTable&,
49  const dimensionSet& dims,
50  const dictionary& dict
51 )
52 :
53  beta_(dict.lookup<scalar>("beta", dimless)),
54  A_(dict.lookup<scalar>("A", dims/pow(dimTemperature, beta_))),
55  Ta_
56  (
57  dict.found("Ta") || !dict.found("Ea")
58  ? dict.lookup<scalar>("Ta", dimTemperature)
59  : dict.lookup<scalar>("Ea", dimEnergy/dimMoles)
60  /constant::physicoChemical::RR.value()
61  ),
62  b_(dict.lookup("b"))
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 {}
70 
71 
73 {}
74 
75 
76 inline Foam::scalar Foam::JanevReactionRate::operator()
77 (
78  const scalar p,
79  const scalar T,
80  const scalarField&,
81  const label
82 ) const
83 {
84  scalar lta = A_;
85 
86  if (mag(beta_) > vSmall)
87  {
88  lta *= pow(T, beta_);
89  }
90 
91  scalar expArg = 0;
92 
93  if (mag(Ta_) > vSmall)
94  {
95  expArg -= Ta_/T;
96  }
97 
98  scalar lnT = log(T);
99 
100  for (int n=0; n<nb_; n++)
101  {
102  expArg += b_[n]*pow(lnT, n);
103  }
104 
105  lta *= exp(expArg);
106 
107  return lta;
108 }
109 
110 
111 inline Foam::scalar Foam::JanevReactionRate::ddT
112 (
113  const scalar p,
114  const scalar T,
115  const scalarField&,
116  const label
117 ) const
118 {
119  scalar lta = A_;
120 
121  if (mag(beta_) > vSmall)
122  {
123  lta *= pow(T, beta_);
124  }
125 
126  scalar expArg = 0;
127 
128  if (mag(Ta_) > vSmall)
129  {
130  expArg -= Ta_/T;
131  }
132 
133  scalar lnT = log(T);
134 
135  for (int n=0; n<nb_; n++)
136  {
137  expArg += b_[n]*pow(lnT, n);
138  }
139 
140  scalar deriv = b_[1];
141 
142  for (int n=2; n<nb_; n++)
143  {
144  deriv += n*b_[n]*pow(lnT, n-1);
145  }
146 
147  lta *= exp(expArg);
148 
149  return lta*(beta_+Ta_/T+deriv)/T;
150 }
151 
152 
154 {
155  return false;
156 }
157 
158 
160 (
161  const scalar p,
162  const scalar T,
163  const scalarField& c,
164  const label li,
165  scalarField& ddc
166 ) const
167 {
168  ddc = 0;
169 }
170 
171 
173 {
174  writeKeyword(os, "A") << A_ << nl;
175  writeKeyword(os, "beta") << beta_ << nl;
176  writeKeyword(os, "Ta") << Ta_ << nl;
177  writeKeyword(os, "b") << b_ << nl;
178 }
179 
180 
181 inline Foam::Ostream& Foam::operator<<
182 (
183  Ostream& os,
184  const JanevReactionRate& jrr
185 )
186 {
187  jrr.write(os);
188  return os;
189 }
190 
191 
192 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
bool found
label n
Janev, Langer, Evans and Post reaction rate.
void preEvaluate() const
Pre-evaluation hook.
void postEvaluate() const
Post-evaluation hook.
void write(Ostream &os) const
Write to stream.
void ddc(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &ddc) const
The derivative of the rate w.r.t. concentration.
JanevReactionRate(const scalar A, const scalar beta, const scalar Ta, const FixedList< scalar, nb_ > b)
Construct from components.
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
The derivative of the rate w.r.t. temperature.
bool hasDdc() const
Is the rate a function of concentration?
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual Ostream & write(const char)=0
Write character.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
A wordList with hashed indices for faster lookup by name.
volScalarField & b
Definition: createFields.H:25
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar exp(const dimensionedScalar &ds)
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
const dimensionSet dimEnergy
const dimensionSet dimless
const dimensionSet dimTemperature
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionSet dimMoles
Ostream & writeKeyword(Foam::Ostream &os, const keyType &kw)
Write the keyword to the Ostream with the current level of indentation.
Definition: keyType.C:155
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p