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-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,
33  const scalar B,
34  const scalar C
35 )
36 :
37  A_(A),
38  beta_(beta),
39  Ta_(Ta),
40  B_(B),
41  C_(C)
42 {}
43 
44 
46 (
47  const speciesTable&,
48  const dictionary& dict
49 )
50 :
51  A_(readScalar(dict.lookup("A"))),
52  beta_(readScalar(dict.lookup("beta"))),
53  Ta_(readScalar(dict.lookup("Ta"))),
54  B_(readScalar(dict.lookup("B"))),
55  C_(readScalar(dict.lookup("C")))
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
60 
61 inline Foam::scalar Foam::LandauTellerReactionRate::operator()
62 (
63  const scalar p,
64  const scalar T,
65  const scalarField&
66 ) const
67 {
68  scalar lta = A_;
69 
70  if (mag(beta_) > vSmall)
71  {
72  lta *= pow(T, beta_);
73  }
74 
75  scalar expArg = 0;
76 
77  if (mag(Ta_) > vSmall)
78  {
79  expArg -= Ta_/T;
80  }
81 
82  if (mag(B_) > vSmall)
83  {
84  expArg += B_/cbrt(T);
85  }
86 
87  if (mag(C_) > vSmall)
88  {
89  expArg += C_/pow(T, 2.0/3.0);
90  }
91 
92  if (mag(expArg) > vSmall)
93  {
94  lta *= exp(expArg);
95  }
96 
97  return lta;
98 }
99 
100 
101 inline Foam::scalar Foam::LandauTellerReactionRate::ddT
102 (
103  const scalar p,
104  const scalar T,
105  const scalarField&
106 ) const
107 {
108  scalar lta = A_;
109 
110  if (mag(beta_) > vSmall)
111  {
112  lta *= pow(T, beta_);
113  }
114 
115  scalar expArg = 0;
116  scalar deriv = 0;
117 
118  if (mag(Ta_) > vSmall)
119  {
120  scalar TaT = Ta_/T;
121  expArg -= TaT;
122  deriv += TaT;
123  }
124 
125  if (mag(B_) > vSmall)
126  {
127  scalar BT = B_/cbrt(T);
128  expArg += BT;
129  deriv -= BT/3;
130  }
131 
132  if (mag(C_) > vSmall)
133  {
134  scalar CT = C_/pow(T, 2.0/3.0);
135  expArg += CT;
136  deriv -= 2*CT/3;
137  }
138 
139  if (mag(expArg) > vSmall)
140  {
141  lta *= exp(expArg);
142  }
143 
144  return lta*(beta_+deriv)/T;
145 }
146 
147 
150 {
151  return NullObjectRef<List<Tuple2<label, scalar>>>();
152 }
153 
154 
156 (
157  const scalar p,
158  const scalar T,
159  const scalarField& c,
160  scalarField& dcidc
161 ) const
162 {}
163 
164 
165 inline Foam::scalar Foam::LandauTellerReactionRate::dcidT
166 (
167  const scalar p,
168  const scalar T,
169  const scalarField& c
170 ) const
171 {
172  return 0;
173 }
174 
175 
177 {
178  writeEntry(os, "A", A_);
179  writeEntry(os, "beta", beta_);
180  writeEntry(os, "Ta", Ta_);
181  writeEntry(os, "B", B_);
182  writeEntry(os, "C", C_);
183 }
184 
185 
186 inline Foam::Ostream& Foam::operator<<
187 (
188  Ostream& os,
189  const LandauTellerReactionRate& arr
190 )
191 {
192  arr.write(os);
193  return os;
194 }
195 
196 
197 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
LandauTellerReactionRate(const scalar A, const scalar beta, const scalar Ta, const scalar B, const scalar C)
Construct from components.
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
void dcidc(const scalar p, const scalar T, const scalarField &c, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
scalar ddT(const scalar p, const scalar T, const scalarField &c) const
const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const volScalarField & T
scalar dcidT(const scalar p, const scalar T, const scalarField &c) const
Temperature derivative of the pressure dependent term.
void write(Ostream &os) const
Write to stream.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Landau-Teller reaction rate.
A wordList with hashed indices for faster lookup by name.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual Ostream & write(const token &)=0
Write next token to stream.
volScalarField & p
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583