TroeFallOffFunctionI.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 "TroeFallOffFunction.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 (
32  const scalar alpha,
33  const scalar Tsss,
34  const scalar Ts,
35  const scalar Tss
36 )
37 :
38  alpha_(alpha),
39  Tsss_(Tsss),
40  Ts_(Ts),
41  Tss_(Tss)
42 {}
43 
44 
46 :
47  alpha_(dict.lookup<scalar>("alpha")),
48  Tsss_(dict.lookup<scalar>("Tsss")),
49  Ts_(dict.lookup<scalar>("Ts")),
50  Tss_(dict.lookup<scalar>("Tss"))
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55 
56 inline Foam::scalar Foam::TroeFallOffFunction::operator()
57 (
58  const scalar T,
59  const scalar Pr
60 ) const
61 {
62  const scalar logPr = log10(max(Pr, small));
63 
64  const scalar Fcent =
65  (1 - alpha_)*exp(-T/Tsss_)
66  + alpha_*exp(-T/Ts_)
67  + exp(-Tss_/T);
68 
69  const scalar logFcent = log10(max(Fcent, small));
70 
71  const scalar c = -0.4 - 0.67*logFcent;
72  static const scalar d = 0.14;
73  const scalar n = 0.75 - 1.27*logFcent;
74 
75  const scalar x1 = n - d*(logPr + c);
76  const scalar x2 = (logPr + c)/x1;
77  const scalar x3 = 1 + sqr(x2);
78  const scalar x4 = logFcent/x3;
79 
80  return pow(10, x4);
81 }
82 
83 
84 inline Foam::scalar Foam::TroeFallOffFunction::ddT
85 (
86  const scalar T,
87  const scalar Pr,
88  const scalar F
89 ) const
90 {
91  static const scalar logTen = log(scalar(10));
92 
93  const scalar logPr = log10(max(Pr, small));
94 
95  const scalar Fcent =
96  (1 - alpha_)*exp(-T/Tsss_)
97  + alpha_*exp(-T/Ts_)
98  + exp(-Tss_/T);
99  const scalar dFcentdT =
100  - (1 - alpha_)/Tsss_*exp(-T/Tsss_)
101  - alpha_/Ts_*exp(-T/Ts_)
102  + Tss_/sqr(T)*exp(-Tss_/T);
103 
104  const scalar logFcent = log10(max(Fcent, small));
105  const scalar dlogFcentdT = Fcent >= small ? dFcentdT/Fcent/logTen : 0;
106 
107  const scalar c = -0.4 - 0.67*logFcent;
108  const scalar dcdT = -0.67*dlogFcentdT;
109  static const scalar d = 0.14;
110  const scalar n = 0.75 - 1.27*logFcent;
111  const scalar dndT = - 1.27*dlogFcentdT;
112 
113  const scalar x1 = n - d*(logPr + c);
114  const scalar dx1dT = dndT - d*dcdT;
115  const scalar x2 = (logPr + c)/x1;
116  const scalar dx2dT = (dcdT - x2*dx1dT)/x1;
117  const scalar x3 = 1 + sqr(x2);
118  const scalar dx3dT = 2*x2*dx2dT;
119  const scalar x4 = logFcent/x3;
120  const scalar dx4dT = (dlogFcentdT - x4*dx3dT)/x3;
121 
122  return logTen*F*dx4dT;
123 }
124 
125 
127 (
128  const scalar T,
129  const scalar Pr,
130  const scalar F
131 ) const
132 {
133  static const scalar logTen = log(scalar(10));
134 
135  const scalar logPr = log10(max(Pr, small));
136  const scalar dlogPrdPr = Pr >= small ? 1/(logTen*Pr) : 0;
137 
138  const scalar Fcent =
139  (1 - alpha_)*exp(-T/Tsss_)
140  + alpha_*exp(-T/Ts_)
141  + exp(-Tss_/T);
142 
143  const scalar logFcent = log10(max(Fcent, small));
144 
145  const scalar c = -0.4 - 0.67*logFcent;
146  static const scalar d = 0.14;
147  const scalar n = 0.75 - 1.27*logFcent;
148 
149  const scalar x1 = n - d*(logPr + c);
150  const scalar dx1dPr = -d*dlogPrdPr;
151  const scalar x2 = (logPr + c)/x1;
152  const scalar dx2dPr = (dlogPrdPr - x2*dx1dPr)/x1;
153  const scalar x3 = 1 + sqr(x2);
154  const scalar dx3dPr = 2*x2*dx2dPr;
155  const scalar x4 = logFcent/x3;
156  const scalar dx4dPr = -x4*dx3dPr/x3;
157 
158  return logTen*F*dx4dPr;
159 }
160 
161 
163 {
164  writeEntry(os, "alpha", alpha_);
165  writeEntry(os, "Tsss", Tsss_);
166  writeEntry(os, "Ts", Ts_);
167  writeEntry(os, "Tss", Tss_);
168 }
169 
170 
171 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
172 
173 inline Foam::Ostream& Foam::operator<<
174 (
175  Foam::Ostream& os,
176  const Foam::TroeFallOffFunction& tfof
177 )
178 {
179  tfof.write(os);
180  return os;
181 }
182 
183 
184 // ************************************************************************* //
label n
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.
The Troe fall-off function.
scalar ddT(const scalar T, const scalar Pr, const scalar F) const
TroeFallOffFunction(const scalar alpha, const scalar Tsss, const scalar Ts, const scalar Tss)
Construct from components.
void write(Ostream &os) const
Write to stream.
scalar ddPr(const scalar T, const scalar Pr, const scalar F) const
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar F
Faraday constant: default SI units: [C/kmol].
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar log10(const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict