LandauTellerReactionRateI.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 
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const scalar A,
34  const scalar beta,
35  const scalar Ta,
36  const scalar B,
37  const scalar C
38 )
39 :
40  beta_(beta),
41  A_(A),
42  Ta_(Ta),
43  B_(B),
44  C_(C)
45 {}
46 
47 
49 (
50  const speciesTable&,
51  const dimensionSet& dims,
52  const dictionary& dict
53 )
54 :
55  beta_(dict.lookup<scalar>("beta", dimless)),
56  A_(dict.lookup<scalar>("A", dims/pow(dimTemperature, beta_))),
57  Ta_
58  (
59  dict.found("Ta") || !dict.found("Ea")
60  ? dict.lookup<scalar>("Ta", dimTemperature)
61  : dict.lookup<scalar>("Ea", dimEnergy/dimMoles)
62  /constant::physicoChemical::RR.value()
63  ),
64  B_(dict.lookup<scalar>("B", cbrt(dimTemperature))),
65  C_(dict.lookup<scalar>("C", sqr(cbrt(dimTemperature))))
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70 
72 {}
73 
74 
76 {}
77 
78 
79 inline Foam::scalar Foam::LandauTellerReactionRate::operator()
80 (
81  const scalar p,
82  const scalar T,
83  const scalarField&,
84  const label
85 ) const
86 {
87  scalar lta = A_;
88 
89  if (mag(beta_) > vSmall)
90  {
91  lta *= pow(T, beta_);
92  }
93 
94  scalar expArg = 0;
95 
96  if (mag(Ta_) > vSmall)
97  {
98  expArg -= Ta_/T;
99  }
100 
101  if (mag(B_) > vSmall)
102  {
103  expArg += B_/cbrt(T);
104  }
105 
106  if (mag(C_) > vSmall)
107  {
108  expArg += C_/sqr(cbrt(T));
109  }
110 
111  if (mag(expArg) > vSmall)
112  {
113  lta *= exp(expArg);
114  }
115 
116  return lta;
117 }
118 
119 
121 (
122  const scalar p,
123  const scalar T,
124  const scalarField&,
125  const label
126 ) const
127 {
128  scalar lta = A_;
129 
130  if (mag(beta_) > vSmall)
131  {
132  lta *= pow(T, beta_);
133  }
134 
135  scalar expArg = 0;
136  scalar deriv = 0;
137 
138  if (mag(Ta_) > vSmall)
139  {
140  scalar TaT = Ta_/T;
141  expArg -= TaT;
142  deriv += TaT;
143  }
144 
145  if (mag(B_) > vSmall)
146  {
147  scalar BT = B_/cbrt(T);
148  expArg += BT;
149  deriv -= BT/3;
150  }
151 
152  if (mag(C_) > vSmall)
153  {
154  scalar CT = C_/sqr(cbrt(T));
155  expArg += CT;
156  deriv -= 2*CT/3;
157  }
158 
159  if (mag(expArg) > vSmall)
160  {
161  lta *= exp(expArg);
162  }
163 
164  return lta*(beta_ + deriv)/T;
165 }
166 
167 
169 {
170  return false;
171 }
172 
173 
175 (
176  const scalar p,
177  const scalar T,
178  const scalarField& c,
179  const label li,
180  scalarField& ddc
181 ) const
182 {
183  ddc = 0;
184 }
185 
186 
188 {
189  writeEntry(os, "A", A_);
190  writeEntry(os, "beta", beta_);
191  writeEntry(os, "Ta", Ta_);
192  writeEntry(os, "B", B_);
193  writeEntry(os, "C", C_);
194 }
195 
196 
197 inline Foam::Ostream& Foam::operator<<
198 (
199  Ostream& os,
200  const LandauTellerReactionRate& arr
201 )
202 {
203  arr.write(os);
204  return os;
205 }
206 
207 
208 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("B", Foam::dimless, 18.678)
bool found
Graphite solid properties.
Definition: C.H:51
Landau-Teller 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.
LandauTellerReactionRate(const scalar A, const scalar beta, const scalar Ta, const scalar B, const scalar C)
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.
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimless
const dimensionSet dimTemperature
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionSet dimMoles
dimensionedScalar cbrt(const dimensionedScalar &ds)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p