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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
27 
28 template<class ReactionRate, class ChemicallyActivationFunction>
30 <
31  ReactionRate,
32  ChemicallyActivationFunction
33 >::ChemicallyActivatedReactionRate
34 (
35  const ReactionRate& k0,
36  const ReactionRate& kInf,
37  const ChemicallyActivationFunction& F,
38  const thirdBodyEfficiencies& tbes
39 )
40 :
41  k0_(k0),
42  kInf_(kInf),
43  F_(F),
44  thirdBodyEfficiencies_(tbes)
45 {
46  forAll(tbes, i)
47  {
48  beta_.append(Tuple2<label, scalar>(i, tbes[i]));
49  }
50 }
51 
52 
53 template<class ReactionRate, class ChemicallyActivationFunction>
55 <
56  ReactionRate,
57  ChemicallyActivationFunction
59 (
60  const speciesTable& species,
61  const dictionary& dict
62 )
63 :
64  k0_(species, dict),
65  kInf_(species, dict),
66  F_(dict),
67  thirdBodyEfficiencies_(species, dict)
68 {
69  forAll(thirdBodyEfficiencies_, i)
70  {
71  beta_.append
72  (
74  (
75  i,
76  thirdBodyEfficiencies_[i]
77  )
78  );
79  }
80 }
81 
82 
83 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
84 
85 template<class ReactionRate, class ChemicallyActivationFunction>
87 <
88  ReactionRate,
89  ChemicallyActivationFunction
90 >::preEvaluate() const
91 {
92  k0_.preEvaluate();
93  kInf_.preEvaluate();
94 }
95 
96 
97 template<class ReactionRate, class ChemicallyActivationFunction>
99 <
100  ReactionRate,
101  ChemicallyActivationFunction
102 >::postEvaluate() const
103 {
104  k0_.postEvaluate();
105  kInf_.postEvaluate();
106 }
107 
108 
109 template<class ReactionRate, class ChemicallyActivationFunction>
110 inline Foam::scalar Foam::ChemicallyActivatedReactionRate
111 <
112  ReactionRate,
113  ChemicallyActivationFunction
114 >::operator()
115 (
116  const scalar p,
117  const scalar T,
118  const scalarField& c,
119  const label li
120 ) const
121 {
122  const scalar k0 = k0_(p, T, c, li);
123  const scalar kInf = kInf_(p, T, c, li);
124  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
125 
126  return k0*(1/(1 + Pr))*F_(T, Pr);
127 }
128 
129 
130 template<class ReactionRate, class ChemicallyActivationFunction>
131 inline Foam::scalar Foam::ChemicallyActivatedReactionRate
132 <
133  ReactionRate,
134  ChemicallyActivationFunction
135 >::ddT
136 (
137  const scalar p,
138  const scalar T,
139  const scalarField& c,
140  const label li
141 ) const
142 {
143  const scalar k0 = k0_(p, T, c, li);
144  const scalar kInf = kInf_(p, T, c, li);
145  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
146 
147  return (1/(1 + Pr))*F_(T, Pr)*k0_.ddT(p, T, c, li);
148 }
149 
150 
151 template<class ReactionRate, class ChemicallyActivationFunction>
153 <
154  ReactionRate,
155  ChemicallyActivationFunction
156 >::dcidc
157 (
158  const scalar p,
159  const scalar T,
160  const scalarField& c,
161  const label li,
162  scalarField& dcidc
163 ) const
164 {
165  const scalar M = thirdBodyEfficiencies_.M(c);
166 
167  if (M > small)
168  {
169  const scalar k0 = k0_(p, T, c, li);
170  const scalar kInf = kInf_(p, T, c, li);
171 
172  const scalar Pr = k0*M/kInf;
173  const scalar F = F_(T, Pr);
174 
175  forAll(beta_, i)
176  {
177  const scalar dPrdci = -beta_[i].second()*k0/kInf;
178  const scalar dFdci = F_.ddc(Pr, F, dPrdci, T);
179  dcidc[i] = (-dPrdci/(1 + Pr) + dFdci/F);
180  }
181  }
182  else
183  {
184  forAll(beta_, i)
185  {
186  dcidc[i] = 0;
187  }
188  }
189 }
190 
191 
192 template<class ReactionRate, class ChemicallyActivationFunction>
193 inline Foam::scalar Foam::ChemicallyActivatedReactionRate
194 <
195  ReactionRate,
196  ChemicallyActivationFunction
197 >::dcidT
198 (
199  const scalar p,
200  const scalar T,
201  const scalarField& c,
202  const label li
203 ) const
204 {
205  const scalar M = thirdBodyEfficiencies_.M(c);
206 
207  if (M > small)
208  {
209  const scalar k0 = k0_(p, T, c, li);
210  const scalar kInf = kInf_(p, T, c, li);
211 
212  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
213  const scalar F = F_(T, Pr);
214  const scalar dPrdT =
215  Pr*(k0_.ddT(p, T, c, li)/k0 - kInf_.ddT(p, T, c, li)/kInf - 1/T);
216  const scalar dFdT = F_.ddT(Pr, F, dPrdT, T);
217 
218  return (-dPrdT/(1 + Pr) + dFdT/F);
219  }
220  else
221  {
222  return 0;
223  }
224 }
225 
226 
227 template<class ReactionRate, class ChemicallyActivationFunction>
229 <
230  ReactionRate,
231  ChemicallyActivationFunction
232 >::write(Ostream& os) const
233 {
234  k0_.write(os);
235  kInf_.write(os);
236  F_.write(os);
237  thirdBodyEfficiencies_.write(os);
238 }
239 
240 
241 template<class ReactionRate, class ChemicallyActivationFunction>
242 inline Foam::Ostream& Foam::operator<<
243 (
244  Ostream& os,
246  <ReactionRate, ChemicallyActivationFunction>& carr
247 )
248 {
249  carr.write(os);
250  return os;
251 }
252 
253 
254 // ************************************************************************* //
dictionary dict
#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
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
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.
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
const volScalarField & T
A wordList with hashed indices for faster lookup by name.
volScalarField & p
Third body efficiencies.
#define M(I)