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-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 // * * * * * * * * * * * * * * * * 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 Foam::scalar
82 (
83  const scalar p,
84  const scalar T,
85  const scalarField& c
86 ) const
87 {
88  const scalar k0 = k0_(p, T, c);
89  const scalar kInf = kInf_(p, T, c);
90  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
91 
92  return kInf*(Pr/(1 + Pr))*F_(T, Pr);
93 }
94 
95 
96 template<class ReactionRate, class FallOffFunction>
97 inline Foam::scalar
99 (
100  const scalar p,
101  const scalar T,
102  const scalarField& c
103 ) const
104 {
105  const scalar k0 = k0_(p, T, c);
106  const scalar kInf = kInf_(p, T, c);
107  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
108 
109  return (Pr/(1 + Pr))*F_(T, Pr)*kInf_.ddT(p, T, c);
110 }
111 
112 
113 template<class ReactionRate, class FallOffFunction>
115 (
116  const scalar p,
117  const scalar T,
118  const scalarField& c,
119  scalarField& dcidc
120 ) const
121 {
122  const scalar M = thirdBodyEfficiencies_.M(c);
123 
124  if (M > small)
125  {
126  const scalar k0 = k0_(p, T, c);
127  const scalar kInf = kInf_(p, T, c);
128  const scalar Pr = k0*M/kInf;
129  const scalar F = F_(T, Pr);
130 
131  forAll(beta_, i)
132  {
133  const scalar dPrdci = -beta_[i].second()*k0/kInf;
134  const scalar dFdci = F_.ddc(Pr, F, dPrdci, T);
135  dcidc[i] = (dPrdci/(Pr*(1 + Pr)) + dFdci/F);
136  }
137  }
138  else
139  {
140  forAll(beta_, i)
141  {
142  dcidc[i] = 0;
143  }
144  }
145 }
146 
147 
148 template<class ReactionRate, class FallOffFunction>
149 inline Foam::scalar
151 (
152  const scalar p,
153  const scalar T,
154  const scalarField& c
155 ) const
156 {
157  const scalar M = thirdBodyEfficiencies_.M(c);
158 
159  if (M > small)
160  {
161  const scalar k0 = k0_(p, T, c);
162  const scalar kInf = kInf_(p, T, c);
163 
164  const scalar Pr = k0*thirdBodyEfficiencies_.M(c)/kInf;
165  const scalar F = F_(T, Pr);
166  const scalar dPrdT =
167  Pr*(k0_.ddT(p, T, c)/k0 - kInf_.ddT(p, T, c)/kInf - 1/T);
168  const scalar dFdT = F_.ddT(Pr, F, dPrdT, T);
169 
170  return (dPrdT/(Pr*(1 + Pr)) + dFdT/F);
171  }
172  else
173  {
174  return 0;
175  }
176 }
177 
178 
179 template<class ReactionRate, class FallOffFunction>
181 (
182  Ostream& os
183 ) const
184 {
185  os << indent << "k0" << nl;
186  os << indent << token::BEGIN_BLOCK << nl;
187  os << incrIndent;
188  k0_.write(os);
189  os << decrIndent;
190  os << indent << token::END_BLOCK << nl;
191 
192  os << indent << "kInf" << nl;
193  os << indent << token::BEGIN_BLOCK << nl;
194  os << incrIndent;
195  kInf_.write(os);
196  os << decrIndent;
197  os << indent << token::END_BLOCK << nl;
198 
199  os << indent << "F" << nl;
200  os << indent << token::BEGIN_BLOCK << nl;
201  os << incrIndent;
202  F_.write(os);
203  os << decrIndent;
204  os << indent << token::END_BLOCK << nl;
205 
206  os << indent << "thirdBodyEfficiencies" << nl;
207  os << indent << token::BEGIN_BLOCK << nl;
208  os << incrIndent;
209  thirdBodyEfficiencies_.write(os);
210  os << decrIndent;
211  os << indent << token::END_BLOCK << nl;
212 }
213 
214 
215 template<class ReactionRate, class FallOffFunction>
216 inline Foam::Ostream& Foam::operator<<
217 (
218  Ostream& os,
220 )
221 {
222  forr.write(os);
223  return os;
224 }
225 
226 
227 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:226
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
General class for handling unimolecular/recombination fall-off reactions.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:699
scalar dcidT(const scalar p, const scalar T, const scalarField &c) const
Temperature derivative of the pressure dependent term.
void write(Ostream &os) const
Write to stream.
volVectorField F(fluid.F())
void dcidc(const scalar p, const scalar T, const scalarField &c, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:265
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:240
scalar ddT(const scalar p, const scalar T, const scalarField &c) const
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 > &)
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.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:233
#define M(I)