FallOffReactionRateI.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 "FallOffReactionRate.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 #define FALL_OFF_FUNCTION_JACOBIAN false
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class ReactionRate, class FallOffFunction>
37 (
38  const ReactionRate& k0,
39  const ReactionRate& kInf,
40  const FallOffFunction& F,
41  const thirdBodyEfficiencies& tbes
42 )
43 :
44  k0_(k0),
45  kInf_(kInf),
46  F_(F),
47  thirdBodyEfficiencies_(tbes)
48 {}
49 
50 
51 template<class ReactionRate, class FallOffFunction>
54 (
55  const speciesTable& species,
56  const dimensionSet& dims,
57  const dictionary& dict
58 )
59 :
60  k0_(species, dims, dict.subDict("k0")),
61  kInf_(species, dims, dict.subDict("kInf")),
62  F_(dict.subDict("F")),
63  thirdBodyEfficiencies_(species, dict.subDict("thirdBodyEfficiencies"))
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
69 template<class ReactionRate, class FallOffFunction>
70 inline void
72 {
73  k0_.preEvaluate();
74  kInf_.preEvaluate();
75 }
76 
77 
78 template<class ReactionRate, class FallOffFunction>
79 inline void
81 {
82  k0_.postEvaluate();
83  kInf_.postEvaluate();
84 }
85 
86 
87 template<class ReactionRate, class FallOffFunction>
88 inline Foam::scalar
90 (
91  const scalar p,
92  const scalar T,
93  const scalarField& c,
94  const label li
95 ) const
96 {
97  const scalar k0 = k0_(p, T, c, li);
98  const scalar kInf = kInf_(p, T, c, li);
99  const scalar M = thirdBodyEfficiencies_.M(c);
100  const scalar Pr = k0/kInf*M;
101  const scalar F = F_(T, Pr);
102 
103  return kInf*(Pr/(1 + Pr))*F;
104 }
105 
106 
107 template<class ReactionRate, class FallOffFunction>
108 inline Foam::scalar
110 (
111  const scalar p,
112  const scalar T,
113  const scalarField& c,
114  const label li
115 ) const
116 {
117  const scalar k0 = k0_(p, T, c, li);
118  const scalar kInf = kInf_(p, T, c, li);
119  const scalar M = thirdBodyEfficiencies_.M(c);
120  const scalar Pr = k0/kInf*M;
121  const scalar F = F_(T, Pr);
122 
123  const scalar dkInfdT = kInf_.ddT(p, T, c, li);
124 
125  #if FALL_OFF_FUNCTION_JACOBIAN
126 
127  const scalar dk0dT = k0_.ddT(p, T, c, li);
128  const scalar dPrdT = (M*dk0dT - Pr*dkInfdT)/kInf;
129  const scalar dFdT = F_.ddT(T, Pr, F) + F_.ddPr(T, Pr, F)*dPrdT;
130 
131  return
132  dkInfdT*(Pr/(1 + Pr))*F
133  + kInf*dPrdT/sqr(1 + Pr)*F
134  + kInf*(Pr/(1 + Pr))*dFdT;
135 
136  #else
137 
138  return dkInfdT*(Pr/(1 + Pr))*F;
139 
140  #endif
141 }
142 
143 
144 template<class ReactionRate, class FallOffFunction>
145 inline bool
147 {
148  return true;
149 }
150 
151 
152 template<class ReactionRate, class FallOffFunction>
154 (
155  const scalar p,
156  const scalar T,
157  const scalarField& c,
158  const label li,
159  scalarField& ddc
160 ) const
161 {
162  const scalar k0 = k0_(p, T, c, li);
163  const scalar kInf = kInf_(p, T, c, li);
164  const scalar M = thirdBodyEfficiencies_.M(c);
165  const scalar Pr = k0/kInf*M;
166  const scalar F = F_(T, Pr);
167 
168  #if FALL_OFF_FUNCTION_JACOBIAN
169 
170  scalarField dk0dc(c.size(), 0);
171  k0_.ddc(p, T, c, li, dk0dc);
172  scalarField dkInfdc(c.size(), 0);
173  kInf_.ddc(p, T, c, li, dkInfdc);
174  tmp<scalarField> tdMdc = thirdBodyEfficiencies_.dMdc(c);
175  const scalarField& dMdc = tdMdc();
176  scalarField dPrdc((M*dk0dc - Pr*dkInfdc + k0*dMdc)/kInf);
177  const scalar dFdPr = F_.ddPr(T, Pr, F);
178 
179  ddc =
180  dkInfdc*(Pr/(1 + Pr))*F
181  + kInf*dPrdc/sqr(1 + Pr)*F
182  + kInf*(Pr/(1 + Pr))*dFdPr*dPrdc;
183 
184  #else
185 
186  kInf_.ddc(p, T, c, li, ddc);
187 
188  ddc *= (Pr/(1 + Pr))*F;
189 
190  #endif
191 }
192 
193 
194 template<class ReactionRate, class FallOffFunction>
196 (
197  Ostream& os
198 ) const
199 {
200  os << indent << "k0" << nl;
201  os << indent << token::BEGIN_BLOCK << nl;
202  os << incrIndent;
203  k0_.write(os);
204  os << decrIndent;
205  os << indent << token::END_BLOCK << nl;
206 
207  os << indent << "kInf" << nl;
208  os << indent << token::BEGIN_BLOCK << nl;
209  os << incrIndent;
210  kInf_.write(os);
211  os << decrIndent;
212  os << indent << token::END_BLOCK << nl;
213 
214  os << indent << "F" << nl;
215  os << indent << token::BEGIN_BLOCK << nl;
216  os << incrIndent;
217  F_.write(os);
218  os << decrIndent;
219  os << indent << token::END_BLOCK << nl;
220 
221  os << indent << "thirdBodyEfficiencies" << nl;
222  os << indent << token::BEGIN_BLOCK << nl;
223  os << incrIndent;
224  thirdBodyEfficiencies_.write(os);
225  os << decrIndent;
226  os << indent << token::END_BLOCK << nl;
227 }
228 
229 
230 template<class ReactionRate, class FallOffFunction>
231 inline Foam::Ostream& Foam::operator<<
232 (
233  Ostream& os,
235 )
236 {
237  forr.write(os);
238  return os;
239 }
240 
241 
242 // ************************************************************************* //
#define M(I)
General class for handling unimolecular/recombination fall-off reactions.
void preEvaluate() const
Pre-evaluation hook.
FallOffReactionRate(const ReactionRate &k0, const ReactionRate &kInf, const FallOffFunction &F, const thirdBodyEfficiencies &tbes)
Construct from components.
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.
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
The derivative of the rate w.r.t. temperature.
friend Ostream & operator(Ostream &, const FallOffReactionRate< ReactionRate, FallOffFunction > &)
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.
Third body efficiencies.
A class for managing temporary objects.
Definition: tmp.H:55
@ BEGIN_BLOCK
Definition: token.H:113
@ END_BLOCK
Definition: token.H:114
const dimensionedScalar F
Faraday constant: default SI units: [C/kmol].
const dimensionedScalar c
Speed of light in a vacuum.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:241
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:234
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p