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-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 #define CHEMICALLY_ACTIVATION_FUNCTION_JACOBIAN false
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 template<class ReactionRate, class ChemicallyActivationFunction>
32 <
33  ReactionRate,
34  ChemicallyActivationFunction
35 >::ChemicallyActivatedReactionRate
36 (
37  const ReactionRate& k0,
38  const ReactionRate& kInf,
39  const ChemicallyActivationFunction& F,
40  const thirdBodyEfficiencies& tbes
41 )
42 :
43  k0_(k0),
44  kInf_(kInf),
45  F_(F),
46  thirdBodyEfficiencies_(tbes)
47 {}
48 
49 
50 template<class ReactionRate, class ChemicallyActivationFunction>
52 <
53  ReactionRate,
54  ChemicallyActivationFunction
56 (
57  const speciesTable& species,
58  const dictionary& dict
59 )
60 :
61  k0_(species, dict.subDict("k0")),
62  kInf_(species, dict.subDict("kInf")),
63  F_(dict.subDict("F")),
64  thirdBodyEfficiencies_(species, dict.subDict("thirdBodyEfficiencies"))
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
70 template<class ReactionRate, class ChemicallyActivationFunction>
72 <
73  ReactionRate,
74  ChemicallyActivationFunction
75 >::preEvaluate() const
76 {
77  k0_.preEvaluate();
78  kInf_.preEvaluate();
79 }
80 
81 
82 template<class ReactionRate, class ChemicallyActivationFunction>
84 <
85  ReactionRate,
86  ChemicallyActivationFunction
87 >::postEvaluate() const
88 {
89  k0_.postEvaluate();
90  kInf_.postEvaluate();
91 }
92 
93 
94 template<class ReactionRate, class ChemicallyActivationFunction>
95 inline Foam::scalar Foam::ChemicallyActivatedReactionRate
96 <
97  ReactionRate,
98  ChemicallyActivationFunction
99 >::operator()
100 (
101  const scalar p,
102  const scalar T,
103  const scalarField& c,
104  const label li
105 ) const
106 {
107  const scalar k0 = k0_(p, T, c, li);
108  const scalar kInf = kInf_(p, T, c, li);
109  const scalar M = thirdBodyEfficiencies_.M(c);
110  const scalar Pr = k0/kInf*M;
111  const scalar F = F_(T, Pr);
112 
113  return k0/(1 + Pr)*F;
114 }
115 
116 
117 template<class ReactionRate, class ChemicallyActivationFunction>
118 inline Foam::scalar Foam::ChemicallyActivatedReactionRate
119 <
120  ReactionRate,
121  ChemicallyActivationFunction
122 >::ddT
123 (
124  const scalar p,
125  const scalar T,
126  const scalarField& c,
127  const label li
128 ) const
129 {
130  const scalar k0 = k0_(p, T, c, li);
131  const scalar kInf = kInf_(p, T, c, li);
132  const scalar M = thirdBodyEfficiencies_.M(c);
133  const scalar Pr = k0/kInf*M;
134  const scalar F = F_(T, Pr);
135 
136  const scalar dk0dT = k0_.ddT(p, T, c, li);
137 
138  #if CHEMICALLY_ACTIVATION_FUNCTION_JACOBIAN
139 
140  const scalar dkInfdT = kInf_.ddT(p, T, c, li);
141  const scalar dPrdT = (M*dk0dT - Pr*dkInfdT)/kInf;
142  const scalar dFdT = F_.ddT(T, Pr, F) + F_.ddPr(T, Pr, F)*dPrdT;
143 
144  return
145  dk0dT/(1 + Pr)*F
146  - k0*dPrdT/sqr(1 + Pr)*F
147  + k0/(1 + Pr)*dFdT;
148 
149  #else
150 
151  return dk0dT/(1 + Pr)*F;
152 
153  #endif
154 }
155 
156 
157 template<class ReactionRate, class ChemicallyActivationFunction>
159 <
160  ReactionRate,
161  ChemicallyActivationFunction
162 >::hasDdc() const
163 {
164  return true;
165 }
166 
167 
168 template<class ReactionRate, class ChemicallyActivationFunction>
170 <
171  ReactionRate,
172  ChemicallyActivationFunction
173 >::ddc
174 (
175  const scalar p,
176  const scalar T,
177  const scalarField& c,
178  const label li,
179  scalarField& ddc
180 ) const
181 {
182  const scalar k0 = k0_(p, T, c, li);
183  const scalar kInf = kInf_(p, T, c, li);
184  const scalar M = thirdBodyEfficiencies_.M(c);
185  const scalar Pr = k0/kInf*M;
186  const scalar F = F_(T, Pr);
187 
188  #if CHEMICALLY_ACTIVATION_FUNCTION_JACOBIAN
189 
190  scalarField dk0dc(c.size(), 0);
191  k0_.ddc(p, T, c, li, dk0dc);
192  scalarField dkInfdc(c.size(), 0);
193  kInf_.ddc(p, T, c, li, dkInfdc);
194  tmp<scalarField> tdMdc = thirdBodyEfficiencies_.dMdc(c);
195  const scalarField& dMdc = tdMdc();
196  scalarField dPrdc((M*dk0dc - Pr*dkInfdc + k0*dMdc)/kInf);
197  const scalar dFdPr = F_.ddPr(T, Pr, F);
198 
199  ddc =
200  dk0dc/(1 + Pr)*F
201  - k0*dPrdc/sqr(1 + Pr)*F
202  + k0/(1 + Pr)*dFdPr*dPrdc;
203 
204  #else
205 
206  k0_.ddc(p, T, c, li, ddc);
207 
208  ddc *= 1/(1 + Pr)*F;
209 
210  #endif
211 }
212 
213 
214 template<class ReactionRate, class ChemicallyActivationFunction>
216 <
217  ReactionRate,
218  ChemicallyActivationFunction
219 >::write(Ostream& os) const
220 {
221  os << indent << "k0" << nl;
222  os << indent << token::BEGIN_BLOCK << nl;
223  os << incrIndent;
224  k0_.write(os);
225  os << decrIndent;
226  os << indent << token::END_BLOCK << nl;
227 
228  os << indent << "kInf" << nl;
229  os << indent << token::BEGIN_BLOCK << nl;
230  os << incrIndent;
231  kInf_.write(os);
232  os << decrIndent;
233  os << indent << token::END_BLOCK << nl;
234 
235  os << indent << "F" << nl;
236  os << indent << token::BEGIN_BLOCK << nl;
237  os << incrIndent;
238  F_.write(os);
239  os << decrIndent;
240  os << indent << token::END_BLOCK << nl;
241 
242  os << indent << "thirdBodyEfficiencies" << nl;
243  os << indent << token::BEGIN_BLOCK << nl;
244  os << incrIndent;
245  thirdBodyEfficiencies_.write(os);
246  os << decrIndent;
247  os << indent << token::END_BLOCK << nl;
248 }
249 
250 
251 template<class ReactionRate, class ChemicallyActivationFunction>
252 inline Foam::Ostream& Foam::operator<<
253 (
254  Ostream& os,
256  <ReactionRate, ChemicallyActivationFunction>& carr
257 )
258 {
259  carr.write(os);
260  return os;
261 }
262 
263 
264 // ************************************************************************* //
virtual Ostream & write(const char)=0
Write character.
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const dimensionedScalar c
Speed of light in a vacuum.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
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:54
static const char nl
Definition: Ostream.H:260
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
const volScalarField & T
A wordList with hashed indices for faster lookup by name.
void postEvaluate() const
Post-evaluation hook.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
Third body efficiencies.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
#define M(I)