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-2021 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 
62 {}
63 
64 
66 {}
67 
68 
69 inline Foam::scalar Foam::LandauTellerReactionRate::operator()
70 (
71  const scalar p,
72  const scalar T,
73  const scalarField&,
74  const label
75 ) const
76 {
77  scalar lta = A_;
78 
79  if (mag(beta_) > vSmall)
80  {
81  lta *= pow(T, beta_);
82  }
83 
84  scalar expArg = 0;
85 
86  if (mag(Ta_) > vSmall)
87  {
88  expArg -= Ta_/T;
89  }
90 
91  if (mag(B_) > vSmall)
92  {
93  expArg += B_/cbrt(T);
94  }
95 
96  if (mag(C_) > vSmall)
97  {
98  expArg += C_/pow(T, 2.0/3.0);
99  }
100 
101  if (mag(expArg) > vSmall)
102  {
103  lta *= exp(expArg);
104  }
105 
106  return lta;
107 }
108 
109 
110 inline Foam::scalar Foam::LandauTellerReactionRate::ddT
111 (
112  const scalar p,
113  const scalar T,
114  const scalarField&,
115  const label
116 ) const
117 {
118  scalar lta = A_;
119 
120  if (mag(beta_) > vSmall)
121  {
122  lta *= pow(T, beta_);
123  }
124 
125  scalar expArg = 0;
126  scalar deriv = 0;
127 
128  if (mag(Ta_) > vSmall)
129  {
130  scalar TaT = Ta_/T;
131  expArg -= TaT;
132  deriv += TaT;
133  }
134 
135  if (mag(B_) > vSmall)
136  {
137  scalar BT = B_/cbrt(T);
138  expArg += BT;
139  deriv -= BT/3;
140  }
141 
142  if (mag(C_) > vSmall)
143  {
144  scalar CT = C_/pow(T, 2.0/3.0);
145  expArg += CT;
146  deriv -= 2*CT/3;
147  }
148 
149  if (mag(expArg) > vSmall)
150  {
151  lta *= exp(expArg);
152  }
153 
154  return lta*(beta_+deriv)/T;
155 }
156 
157 
159 {
160  return false;
161 }
162 
163 
165 (
166  const scalar p,
167  const scalar T,
168  const scalarField& c,
169  const label li,
170  scalarField& ddc
171 ) const
172 {
173  ddc = 0;
174 }
175 
176 
178 {
179  writeEntry(os, "A", A_);
180  writeEntry(os, "beta", beta_);
181  writeEntry(os, "Ta", Ta_);
182  writeEntry(os, "B", B_);
183  writeEntry(os, "C", C_);
184 }
185 
186 
187 inline Foam::Ostream& Foam::operator<<
188 (
189  Ostream& os,
190  const LandauTellerReactionRate& arr
191 )
192 {
193  arr.write(os);
194  return os;
195 }
196 
197 
198 // ************************************************************************* //
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.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
bool hasDdc() const
Is the rate a function of concentration?
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
void postEvaluate() const
Post-evaluation hook.
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.
void preEvaluate() const
Pre-evaluation hook.
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 T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
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
The derivative of the rate w.r.t. temperature.
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:864