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_(dict.lookup<scalar>("A")),
52  beta_(dict.lookup<scalar>("beta")),
53  Ta_(dict.lookup<scalar>("Ta")),
54  B_(dict.lookup<scalar>("B")),
55  C_(dict.lookup<scalar>("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 label
67 ) const
68 {
69  scalar lta = A_;
70 
71  if (mag(beta_) > vSmall)
72  {
73  lta *= pow(T, beta_);
74  }
75 
76  scalar expArg = 0;
77 
78  if (mag(Ta_) > vSmall)
79  {
80  expArg -= Ta_/T;
81  }
82 
83  if (mag(B_) > vSmall)
84  {
85  expArg += B_/cbrt(T);
86  }
87 
88  if (mag(C_) > vSmall)
89  {
90  expArg += C_/pow(T, 2.0/3.0);
91  }
92 
93  if (mag(expArg) > vSmall)
94  {
95  lta *= exp(expArg);
96  }
97 
98  return lta;
99 }
100 
101 
102 inline Foam::scalar Foam::LandauTellerReactionRate::ddT
103 (
104  const scalar p,
105  const scalar T,
106  const scalarField&,
107  const label
108 ) const
109 {
110  scalar lta = A_;
111 
112  if (mag(beta_) > vSmall)
113  {
114  lta *= pow(T, beta_);
115  }
116 
117  scalar expArg = 0;
118  scalar deriv = 0;
119 
120  if (mag(Ta_) > vSmall)
121  {
122  scalar TaT = Ta_/T;
123  expArg -= TaT;
124  deriv += TaT;
125  }
126 
127  if (mag(B_) > vSmall)
128  {
129  scalar BT = B_/cbrt(T);
130  expArg += BT;
131  deriv -= BT/3;
132  }
133 
134  if (mag(C_) > vSmall)
135  {
136  scalar CT = C_/pow(T, 2.0/3.0);
137  expArg += CT;
138  deriv -= 2*CT/3;
139  }
140 
141  if (mag(expArg) > vSmall)
142  {
143  lta *= exp(expArg);
144  }
145 
146  return lta*(beta_+deriv)/T;
147 }
148 
149 
152 {
153  return NullObjectRef<List<Tuple2<label, scalar>>>();
154 }
155 
156 
158 (
159  const scalar p,
160  const scalar T,
161  const scalarField& c,
162  const label li,
163  scalarField& dcidc
164 ) const
165 {}
166 
167 
168 inline Foam::scalar Foam::LandauTellerReactionRate::dcidT
169 (
170  const scalar p,
171  const scalar T,
172  const scalarField& c,
173  const label li
174 ) const
175 {
176  return 0;
177 }
178 
179 
181 {
182  writeEntry(os, "A", A_);
183  writeEntry(os, "beta", beta_);
184  writeEntry(os, "Ta", Ta_);
185  writeEntry(os, "B", B_);
186  writeEntry(os, "C", C_);
187 }
188 
189 
190 inline Foam::Ostream& Foam::operator<<
191 (
192  Ostream& os,
193  const LandauTellerReactionRate& arr
194 )
195 {
196  arr.write(os);
197  return os;
198 }
199 
200 
201 // ************************************************************************* //
dictionary dict
LandauTellerReactionRate(const scalar A, const scalar beta, const scalar Ta, const scalar B, const scalar C)
Construct from components.
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
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, const label li, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
scalar dcidT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of the pressure dependent term.
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:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const volScalarField & T
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.
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812