ChemicallyActivatedReactionRateI.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 
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 #define CHEMICALLY_ACTIVATION_FUNCTION_JACOBIAN false
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class ReactionRate, class ChemicallyActivationFunction>
36 <
37  ReactionRate,
38  ChemicallyActivationFunction
40 (
41  const ReactionRate& k0,
42  const ReactionRate& kInf,
43  const ChemicallyActivationFunction& F,
44  const thirdBodyEfficiencies& tbes
45 )
46 :
47  k0_(k0),
48  kInf_(kInf),
49  F_(F),
50  thirdBodyEfficiencies_(tbes)
51 {}
52 
53 
54 template<class ReactionRate, class ChemicallyActivationFunction>
56 <
57  ReactionRate,
58  ChemicallyActivationFunction
60 (
61  const speciesTable& species,
62  const dimensionSet& dims,
63  const dictionary& dict
64 )
65 :
66  k0_(species, dims, dict.subDict("k0")),
67  kInf_(species, dims, dict.subDict("kInf")),
68  F_(dict.subDict("F")),
69  thirdBodyEfficiencies_(species, dict.subDict("thirdBodyEfficiencies"))
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 template<class ReactionRate, class ChemicallyActivationFunction>
77 <
78  ReactionRate,
79  ChemicallyActivationFunction
80 >::preEvaluate() const
81 {
82  k0_.preEvaluate();
83  kInf_.preEvaluate();
84 }
85 
86 
87 template<class ReactionRate, class ChemicallyActivationFunction>
89 <
90  ReactionRate,
91  ChemicallyActivationFunction
92 >::postEvaluate() const
93 {
94  k0_.postEvaluate();
95  kInf_.postEvaluate();
96 }
97 
98 
99 template<class ReactionRate, class ChemicallyActivationFunction>
100 inline Foam::scalar Foam::ChemicallyActivatedReactionRate
101 <
102  ReactionRate,
103  ChemicallyActivationFunction
104 >::operator()
105 (
106  const scalar p,
107  const scalar T,
108  const scalarField& c,
109  const label li
110 ) const
111 {
112  const scalar k0 = k0_(p, T, c, li);
113  const scalar kInf = kInf_(p, T, c, li);
114  const scalar M = thirdBodyEfficiencies_.M(c);
115  const scalar Pr = k0/kInf*M;
116  const scalar F = F_(T, Pr);
117 
118  return k0/(1 + Pr)*F;
119 }
120 
121 
122 template<class ReactionRate, class ChemicallyActivationFunction>
123 inline Foam::scalar Foam::ChemicallyActivatedReactionRate
124 <
125  ReactionRate,
126  ChemicallyActivationFunction
127 >::ddT
128 (
129  const scalar p,
130  const scalar T,
131  const scalarField& c,
132  const label li
133 ) const
134 {
135  const scalar k0 = k0_(p, T, c, li);
136  const scalar kInf = kInf_(p, T, c, li);
137  const scalar M = thirdBodyEfficiencies_.M(c);
138  const scalar Pr = k0/kInf*M;
139  const scalar F = F_(T, Pr);
140 
141  const scalar dk0dT = k0_.ddT(p, T, c, li);
142 
143  #if CHEMICALLY_ACTIVATION_FUNCTION_JACOBIAN
144 
145  const scalar dkInfdT = kInf_.ddT(p, T, c, li);
146  const scalar dPrdT = (M*dk0dT - Pr*dkInfdT)/kInf;
147  const scalar dFdT = F_.ddT(T, Pr, F) + F_.ddPr(T, Pr, F)*dPrdT;
148 
149  return
150  dk0dT/(1 + Pr)*F
151  - k0*dPrdT/sqr(1 + Pr)*F
152  + k0/(1 + Pr)*dFdT;
153 
154  #else
155 
156  return dk0dT/(1 + Pr)*F;
157 
158  #endif
159 }
160 
161 
162 template<class ReactionRate, class ChemicallyActivationFunction>
164 <
165  ReactionRate,
166  ChemicallyActivationFunction
167 >::hasDdc() const
168 {
169  return true;
170 }
171 
172 
173 template<class ReactionRate, class ChemicallyActivationFunction>
175 <
176  ReactionRate,
177  ChemicallyActivationFunction
178 >::ddc
179 (
180  const scalar p,
181  const scalar T,
182  const scalarField& c,
183  const label li,
184  scalarField& ddc
185 ) const
186 {
187  const scalar k0 = k0_(p, T, c, li);
188  const scalar kInf = kInf_(p, T, c, li);
189  const scalar M = thirdBodyEfficiencies_.M(c);
190  const scalar Pr = k0/kInf*M;
191  const scalar F = F_(T, Pr);
192 
193  #if CHEMICALLY_ACTIVATION_FUNCTION_JACOBIAN
194 
195  scalarField dk0dc(c.size(), 0);
196  k0_.ddc(p, T, c, li, dk0dc);
197  scalarField dkInfdc(c.size(), 0);
198  kInf_.ddc(p, T, c, li, dkInfdc);
199  tmp<scalarField> tdMdc = thirdBodyEfficiencies_.dMdc(c);
200  const scalarField& dMdc = tdMdc();
201  scalarField dPrdc((M*dk0dc - Pr*dkInfdc + k0*dMdc)/kInf);
202  const scalar dFdPr = F_.ddPr(T, Pr, F);
203 
204  ddc =
205  dk0dc/(1 + Pr)*F
206  - k0*dPrdc/sqr(1 + Pr)*F
207  + k0/(1 + Pr)*dFdPr*dPrdc;
208 
209  #else
210 
211  k0_.ddc(p, T, c, li, ddc);
212 
213  ddc *= 1/(1 + Pr)*F;
214 
215  #endif
216 }
217 
218 
219 template<class ReactionRate, class ChemicallyActivationFunction>
221 <
222  ReactionRate,
223  ChemicallyActivationFunction
224 >::write(Ostream& os) const
225 {
226  os << indent << "k0" << nl;
227  os << indent << token::BEGIN_BLOCK << nl;
228  os << incrIndent;
229  k0_.write(os);
230  os << decrIndent;
231  os << indent << token::END_BLOCK << nl;
232 
233  os << indent << "kInf" << nl;
234  os << indent << token::BEGIN_BLOCK << nl;
235  os << incrIndent;
236  kInf_.write(os);
237  os << decrIndent;
238  os << indent << token::END_BLOCK << nl;
239 
240  os << indent << "F" << nl;
241  os << indent << token::BEGIN_BLOCK << nl;
242  os << incrIndent;
243  F_.write(os);
244  os << decrIndent;
245  os << indent << token::END_BLOCK << nl;
246 
247  os << indent << "thirdBodyEfficiencies" << nl;
248  os << indent << token::BEGIN_BLOCK << nl;
249  os << incrIndent;
250  thirdBodyEfficiencies_.write(os);
251  os << decrIndent;
252  os << indent << token::END_BLOCK << nl;
253 }
254 
255 
256 template<class ReactionRate, class ChemicallyActivationFunction>
257 inline Foam::Ostream& Foam::operator<<
258 (
259  Ostream& os,
261  <ReactionRate, ChemicallyActivationFunction>& carr
262 )
263 {
264  carr.write(os);
265  return os;
266 }
267 
268 
269 // ************************************************************************* //
#define M(I)
General class for handling chemically-activated bimolecular reactions.
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.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
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