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-2019 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 label li
96 ) const
97 {
98  const scalar k0 = k0_(p, T, c, li);
99  const scalar kInf = kInf_(p, T, c, li);
100  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
101 
102  return k0*(1/(1 + Pr))*F_(T, Pr);
103 }
104 
105 
106 template<class ReactionRate, class ChemicallyActivationFunction>
107 inline Foam::scalar Foam::ChemicallyActivatedReactionRate
108 <
109  ReactionRate,
110  ChemicallyActivationFunction
111 >::ddT
112 (
113  const scalar p,
114  const scalar T,
115  const scalarField& c,
116  const label li
117 ) const
118 {
119  const scalar k0 = k0_(p, T, c, li);
120  const scalar kInf = kInf_(p, T, c, li);
121  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
122 
123  return (1/(1 + Pr))*F_(T, Pr)*k0_.ddT(p, T, c, li);
124 }
125 
126 
127 template<class ReactionRate, class ChemicallyActivationFunction>
129 <
130  ReactionRate,
131  ChemicallyActivationFunction
132 >::dcidc
133 (
134  const scalar p,
135  const scalar T,
136  const scalarField& c,
137  const label li,
138  scalarField& dcidc
139 ) const
140 {
141  const scalar M = thirdBodyEfficiencies_.M(c);
142 
143  if (M > small)
144  {
145  const scalar k0 = k0_(p, T, c, li);
146  const scalar kInf = kInf_(p, T, c, li);
147 
148  const scalar Pr = k0*M/kInf;
149  const scalar F = F_(T, Pr);
150 
151  forAll(beta_, i)
152  {
153  const scalar dPrdci = -beta_[i].second()*k0/kInf;
154  const scalar dFdci = F_.ddc(Pr, F, dPrdci, T);
155  dcidc[i] = (-dPrdci/(1 + Pr) + dFdci/F);
156  }
157  }
158  else
159  {
160  forAll(beta_, i)
161  {
162  dcidc[i] = 0;
163  }
164  }
165 }
166 
167 
168 template<class ReactionRate, class ChemicallyActivationFunction>
169 inline Foam::scalar Foam::ChemicallyActivatedReactionRate
170 <
171  ReactionRate,
172  ChemicallyActivationFunction
173 >::dcidT
174 (
175  const scalar p,
176  const scalar T,
177  const scalarField& c,
178  const label li
179 ) const
180 {
181  const scalar M = thirdBodyEfficiencies_.M(c);
182 
183  if (M > small)
184  {
185  const scalar k0 = k0_(p, T, c, li);
186  const scalar kInf = kInf_(p, T, c, li);
187 
188  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
189  const scalar F = F_(T, Pr);
190  const scalar dPrdT =
191  Pr*(k0_.ddT(p, T, c, li)/k0 - kInf_.ddT(p, T, c, li)/kInf - 1/T);
192  const scalar dFdT = F_.ddT(Pr, F, dPrdT, T);
193 
194  return (-dPrdT/(1 + Pr) + dFdT/F);
195  }
196  else
197  {
198  return 0;
199  }
200 }
201 
202 
203 template<class ReactionRate, class ChemicallyActivationFunction>
205 <
206  ReactionRate,
207  ChemicallyActivationFunction
208 >::write(Ostream& os) const
209 {
210  k0_.write(os);
211  kInf_.write(os);
212  F_.write(os);
213  thirdBodyEfficiencies_.write(os);
214 }
215 
216 
217 template<class ReactionRate, class ChemicallyActivationFunction>
218 inline Foam::Ostream& Foam::operator<<
219 (
220  Ostream& os,
222  <ReactionRate, ChemicallyActivationFunction>& carr
223 )
224 {
225  carr.write(os);
226  return os;
227 }
228 
229 
230 // ************************************************************************* //
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:158
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)