ReversibleReaction.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 #include "ReversibleReaction.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template
31 <
32  template<class> class ReactionType,
33  class ReactionThermo,
34  class ReactionRate
35 >
38 (
39  const ReactionType<ReactionThermo>& reaction,
40  const ReactionRate& k
41 )
42 :
43  ReactionType<ReactionThermo>(reaction),
44  k_(k)
45 {}
46 
47 
48 template
49 <
50  template<class> class ReactionType,
51  class ReactionThermo,
52  class ReactionRate
53 >
56 (
57  const speciesTable& species,
58  const HashPtrTable<ReactionThermo>& thermoDatabase,
59  Istream& is
60 )
61 :
62  ReactionType<ReactionThermo>(species, thermoDatabase, is),
63  k_(species, is)
64 {}
65 
66 
67 template
68 <
69  template<class> class ReactionType,
70  class ReactionThermo,
71  class ReactionRate
72 >
75 (
76  const speciesTable& species,
77  const HashPtrTable<ReactionThermo>& thermoDatabase,
78  const dictionary& dict
79 )
80 :
81  ReactionType<ReactionThermo>(species, thermoDatabase, dict),
82  k_(species, dict)
83 {}
84 
85 
86 template
87 <
88  template<class> class ReactionType,
89  class ReactionThermo,
90  class ReactionRate
91 >
94 (
96  const speciesTable& species
97 )
98 :
99  ReactionType<ReactionThermo>(rr, species),
100  k_(rr.k_)
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
106 template
107 <
108  template<class> class ReactionType,
109  class ReactionThermo,
110  class ReactionRate
111 >
112 Foam::scalar Foam::ReversibleReaction
113 <
114  ReactionType,
115  ReactionThermo,
116  ReactionRate
117 >::kf
118 (
119  const scalar p,
120  const scalar T,
121  const scalarField& c
122 ) const
123 {
124  return k_(p, T, c);
125 }
126 
127 
128 template
129 <
130  template<class> class ReactionType,
131  class ReactionThermo,
132  class ReactionRate
133 >
134 Foam::scalar Foam::ReversibleReaction
135 <
136  ReactionType,
137  ReactionThermo,
138  ReactionRate
139 >::kr
140 (
141  const scalar kfwd,
142  const scalar p,
143  const scalar T,
144  const scalarField& c
145 ) const
146 {
147  scalar Kc = this->Kc(p, T);
148 
149  if (mag(Kc) > VSMALL)
150  {
151  return kfwd/Kc;
152  }
153  else
154  {
155  return 0;
156  }
157 }
158 
159 
160 template
161 <
162  template<class> class ReactionType,
163  class ReactionThermo,
164  class ReactionRate
165 >
166 Foam::scalar Foam::ReversibleReaction
167 <
168  ReactionType,
169  ReactionThermo,
170  ReactionRate
171 >::kr
172 (
173  const scalar p,
174  const scalar T,
175  const scalarField& c
176 ) const
177 {
178  return kr(kf(p, T, c), p, T, c);
179 }
180 
181 
182 template
183 <
184  template<class> class ReactionType,
185  class ReactionThermo,
186  class ReactionRate
187 >
189 <
190  ReactionType,
191  ReactionThermo,
192  ReactionRate
193 >::write
194 (
195  Ostream& os
196 ) const
197 {
199  k_.write(os);
200 }
201 
202 
203 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
ReversibleReaction(const ReactionType< ReactionThermo > &reaction, const ReactionRate &k)
Construct from components.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A wordList with hashed indices for faster lookup by name.
Info<< "Creating reaction model\n"<< endl;autoPtr< combustionModels::psiCombustionModel > reaction(combustionModels::psiCombustionModel::New(mesh))
dimensioned< scalar > mag(const dimensioned< Type > &)
runTime write()
Simple extension of Reaction to handle reversible reactions using equilibrium thermodynamics.