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-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 
28 template<class ReactionRate, class FallOffFunction>
31 (
32  const ReactionRate& k0,
33  const ReactionRate& kInf,
34  const FallOffFunction& F,
35  const thirdBodyEfficiencies& tbes
36 )
37 :
38  k0_(k0),
39  kInf_(kInf),
40  F_(F),
41  thirdBodyEfficiencies_(tbes)
42 {
43  forAll(tbes, i)
44  {
45  beta_.append(Tuple2<label, scalar>(i, tbes[i]));
46  }
47 }
48 
49 
50 template<class ReactionRate, class FallOffFunction>
53 (
54  const speciesTable& species,
55  const dictionary& dict
56 )
57 :
58  k0_(species, dict.subDict("k0")),
59  kInf_(species, dict.subDict("kInf")),
60  F_(dict.subDict("F")),
61  thirdBodyEfficiencies_(species, dict.subDict("thirdBodyEfficiencies"))
62 {
63  forAll(thirdBodyEfficiencies_, i)
64  {
65  beta_.append
66  (
68  (
69  i,
70  thirdBodyEfficiencies_[i]
71  )
72  );
73  }
74 }
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
79 template<class ReactionRate, class FallOffFunction>
80 inline void
82 {
83  k0_.preEvaluate();
84  kInf_.preEvaluate();
85 }
86 
87 
88 template<class ReactionRate, class FallOffFunction>
89 inline void
91 {
92  k0_.postEvaluate();
93  kInf_.postEvaluate();
94 }
95 
96 
97 template<class ReactionRate, class FallOffFunction>
98 inline Foam::scalar
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 Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
110 
111  return kInf*(Pr/(1 + Pr))*F_(T, Pr);
112 }
113 
114 
115 template<class ReactionRate, class FallOffFunction>
116 inline Foam::scalar
118 (
119  const scalar p,
120  const scalar T,
121  const scalarField& c,
122  const label li
123 ) const
124 {
125  const scalar k0 = k0_(p, T, c, li);
126  const scalar kInf = kInf_(p, T, c, li);
127  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
128 
129  return (Pr/(1 + Pr))*F_(T, Pr)*kInf_.ddT(p, T, c, li);
130 }
131 
132 
133 template<class ReactionRate, class FallOffFunction>
135 (
136  const scalar p,
137  const scalar T,
138  const scalarField& c,
139  const label li,
140  scalarField& dcidc
141 ) const
142 {
143  const scalar M = thirdBodyEfficiencies_.M(c);
144 
145  if (M > small)
146  {
147  const scalar k0 = k0_(p, T, c, li);
148  const scalar kInf = kInf_(p, T, c, li);
149  const scalar Pr = k0*M/kInf;
150  const scalar F = F_(T, Pr);
151 
152  forAll(beta_, i)
153  {
154  const scalar dPrdci = -beta_[i].second()*k0/kInf;
155  const scalar dFdci = F_.ddc(Pr, F, dPrdci, T);
156  dcidc[i] = (dPrdci/(Pr*(1 + Pr)) + dFdci/F);
157  }
158  }
159  else
160  {
161  forAll(beta_, i)
162  {
163  dcidc[i] = 0;
164  }
165  }
166 }
167 
168 
169 template<class ReactionRate, class FallOffFunction>
170 inline Foam::scalar
172 (
173  const scalar p,
174  const scalar T,
175  const scalarField& c,
176  const label li
177 ) const
178 {
179  const scalar M = thirdBodyEfficiencies_.M(c);
180 
181  if (M > small)
182  {
183  const scalar k0 = k0_(p, T, c, li);
184  const scalar kInf = kInf_(p, T, c, li);
185 
186  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
187  const scalar F = F_(T, Pr);
188  const scalar dPrdT =
189  Pr*(k0_.ddT(p, T, c, li)/k0 - kInf_.ddT(p, T, c, li)/kInf - 1/T);
190  const scalar dFdT = F_.ddT(Pr, F, dPrdT, T);
191 
192  return (dPrdT/(Pr*(1 + Pr)) + dFdT/F);
193  }
194  else
195  {
196  return 0;
197  }
198 }
199 
200 
201 template<class ReactionRate, class FallOffFunction>
203 (
204  Ostream& os
205 ) const
206 {
207  os << indent << "k0" << nl;
208  os << indent << token::BEGIN_BLOCK << nl;
209  os << incrIndent;
210  k0_.write(os);
211  os << decrIndent;
212  os << indent << token::END_BLOCK << nl;
213 
214  os << indent << "kInf" << nl;
215  os << indent << token::BEGIN_BLOCK << nl;
216  os << incrIndent;
217  kInf_.write(os);
218  os << decrIndent;
219  os << indent << token::END_BLOCK << nl;
220 
221  os << indent << "F" << nl;
222  os << indent << token::BEGIN_BLOCK << nl;
223  os << incrIndent;
224  F_.write(os);
225  os << decrIndent;
226  os << indent << token::END_BLOCK << nl;
227 
228  os << indent << "thirdBodyEfficiencies" << nl;
229  os << indent << token::BEGIN_BLOCK << nl;
230  os << incrIndent;
231  thirdBodyEfficiencies_.write(os);
232  os << decrIndent;
233  os << indent << token::END_BLOCK << nl;
234 }
235 
236 
237 template<class ReactionRate, class FallOffFunction>
238 inline Foam::Ostream& Foam::operator<<
239 (
240  Ostream& os,
242 )
243 {
244  forr.write(os);
245  return os;
246 }
247 
248 
249 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
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
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
void postEvaluate() const
Post-evaluation hook.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
A 2-tuple for storing two objects of different types.
Definition: HashTable.H:65
General class for handling unimolecular/recombination fall-off reactions.
const dimensionedScalar c
Speed of light in a vacuum.
void dcidc(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:982
scalar dcidT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of the pressure dependent term.
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
void write(Ostream &os) const
Write to stream.
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
FallOffReactionRate(const ReactionRate &k0, const ReactionRate &kInf, const FallOffFunction &F, const thirdBodyEfficiencies &tbes)
Construct from components.
friend Ostream & operator(Ostream &, const FallOffReactionRate< ReactionRate, FallOffFunction > &)
void preEvaluate() const
Pre-evaluation hook.
A wordList with hashed indices for faster lookup by name.
volScalarField & p
Third body efficiencies.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
#define M(I)