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-2018 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>
86 inline Foam::scalar Foam::ChemicallyActivatedReactionRate
87 <
88  ReactionRate,
89  ChemicallyActivationFunction
90 >::operator()
91 (
92  const scalar p,
93  const scalar T,
94  const scalarField& c
95 ) const
96 {
97  const scalar k0 = k0_(p, T, c);
98  const scalar kInf = kInf_(p, T, c);
99  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
100 
101  return k0*(1/(1 + Pr))*F_(T, Pr);
102 }
103 
104 
105 template<class ReactionRate, class ChemicallyActivationFunction>
106 inline Foam::scalar Foam::ChemicallyActivatedReactionRate
107 <
108  ReactionRate,
109  ChemicallyActivationFunction
110 >::ddT
111 (
112  const scalar p,
113  const scalar T,
114  const scalarField& c
115 ) const
116 {
117  const scalar k0 = k0_(p, T, c);
118  const scalar kInf = kInf_(p, T, c);
119  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
120 
121  return (1/(1 + Pr))*F_(T, Pr)*k0_.ddT(p, T, c);
122 }
123 
124 
125 template<class ReactionRate, class ChemicallyActivationFunction>
127 <
128  ReactionRate,
129  ChemicallyActivationFunction
130 >::dcidc
131 (
132  const scalar p,
133  const scalar T,
134  const scalarField& c,
135  scalarField& dcidc
136 ) const
137 {
138  const scalar M = thirdBodyEfficiencies_.M(c);
139 
140  if (M > small)
141  {
142  const scalar k0 = k0_(p, T, c);
143  const scalar kInf = kInf_(p, T, c);
144 
145  const scalar Pr = k0*M/kInf;
146  const scalar F = F_(T, Pr);
147 
148  forAll(beta_, i)
149  {
150  const scalar dPrdci = -beta_[i].second()*k0/kInf;
151  const scalar dFdci = F_.ddc(Pr, F, dPrdci, T);
152  dcidc[i] = (-dPrdci/(1 + Pr) + dFdci/F);
153  }
154  }
155  else
156  {
157  forAll(beta_, i)
158  {
159  dcidc[i] = 0;
160  }
161  }
162 }
163 
164 
165 template<class ReactionRate, class ChemicallyActivationFunction>
166 inline Foam::scalar Foam::ChemicallyActivatedReactionRate
167 <
168  ReactionRate,
169  ChemicallyActivationFunction
170 >::dcidT
171 (
172  const scalar p,
173  const scalar T,
174  const scalarField& c
175 ) const
176 {
177  const scalar M = thirdBodyEfficiencies_.M(c);
178 
179  if (M > small)
180  {
181  const scalar k0 = k0_(p, T, c);
182  const scalar kInf = kInf_(p, T, c);
183 
184  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
185  const scalar F = F_(T, Pr);
186  const scalar dPrdT =
187  Pr*(k0_.ddT(p, T, c)/k0 - kInf_.ddT(p, T, c)/kInf - 1/T);
188  const scalar dFdT = F_.ddT(Pr, F, dPrdT, T);
189 
190  return (-dPrdT/(1 + Pr) + dFdT/F);
191  }
192  else
193  {
194  return 0;
195  }
196 }
197 
198 
199 template<class ReactionRate, class ChemicallyActivationFunction>
201 <
202  ReactionRate,
203  ChemicallyActivationFunction
204 >::write(Ostream& os) const
205 {
206  k0_.write(os);
207  kInf_.write(os);
208  F_.write(os);
209  thirdBodyEfficiencies_.write(os);
210 }
211 
212 
213 template<class ReactionRate, class ChemicallyActivationFunction>
214 inline Foam::Ostream& Foam::operator<<
215 (
216  Ostream& os,
218  <ReactionRate, ChemicallyActivationFunction>& carr
219 )
220 {
221  carr.write(os);
222  return os;
223 }
224 
225 
226 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
A 2-tuple for storing two objects of different types.
Definition: HashTable.H:65
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
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:53
const volScalarField & T
A wordList with hashed indices for faster lookup by name.
const dimensionedScalar c
Speed of light in a vacuum.
virtual Ostream & write(const token &)=0
Write next token to stream.
volScalarField & p
Third body efficiencies.
#define M(I)