NonEquilibriumReversibleReaction.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-2012 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 
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template
31 <
32  template<class> class ReactionType,
33  class ReactionThermo,
34  class ReactionRate
35 >
37 <
38  ReactionType,
39  ReactionThermo,
40  ReactionRate
41 >::
42 NonEquilibriumReversibleReaction
43 (
44  const ReactionType<ReactionThermo>& reaction,
45  const ReactionRate& forwardReactionRate,
46  const ReactionRate& reverseReactionRate
47 )
48 :
49  ReactionType<ReactionThermo>(reaction),
50  fk_(forwardReactionRate),
51  rk_(reverseReactionRate)
52 {}
53 
54 
55 
56 template
57 <
58  template<class> class ReactionType,
59  class ReactionThermo,
60  class ReactionRate
61 >
63 <
64  ReactionType,
65  ReactionThermo,
66  ReactionRate
67 >::
69 (
70  const speciesTable& species,
71  const HashPtrTable<ReactionThermo>& thermoDatabase,
72  Istream& is
73 )
74 :
75  ReactionType<ReactionThermo>(species, thermoDatabase, is),
76  fk_(species, is),
77  rk_(species, is)
78 {}
79 
80 
81 template
82 <
83  template<class> class ReactionType,
84  class ReactionThermo,
85  class ReactionRate
86 >
88 <
89  ReactionType,
90  ReactionThermo,
91  ReactionRate
92 >::
94 (
95  const speciesTable& species,
96  const HashPtrTable<ReactionThermo>& thermoDatabase,
97  const dictionary& dict
98 )
99 :
100  ReactionType<ReactionThermo>(species, thermoDatabase, dict),
101  fk_(species, dict.subDict("forward")),
102  rk_(species, dict.subDict("reverse"))
103 {}
104 
105 
106 template
107 <
108  template<class> class ReactionType,
109  class ReactionThermo,
110  class ReactionRate
111 >
113 <
114  ReactionType,
115  ReactionThermo,
116  ReactionRate
117 >::
119 (
121  <
122  ReactionType,
123  ReactionThermo,
124  ReactionRate
125  >& nerr,
126  const speciesTable& species
127 )
128 :
129  ReactionType<ReactionThermo>(nerr, species),
130  fk_(nerr.fk_),
131  rk_(nerr.rk_)
132 {}
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
137 template
138 <
139  template<class> class ReactionType,
140  class ReactionThermo,
141  class ReactionRate
142 >
143 Foam::scalar
145 <
146  ReactionType,
147  ReactionThermo,
148  ReactionRate
149 >::kf
150 (
151  const scalar p,
152  const scalar T,
153  const scalarField& c
154 ) const
155 {
156  return fk_(p, T, c);
157 }
158 
159 
160 template
161 <
162  template<class> class ReactionType,
163  class ReactionThermo,
164  class ReactionRate
165 >
166 Foam::scalar
168 <
169  ReactionType,
170  ReactionThermo,
171  ReactionRate
172 >::kr
173 (
174  const scalar,
175  const scalar p,
176  const scalar T,
177  const scalarField& c
178 ) const
179 {
180  return rk_(p, T, c);
181 }
182 
183 
184 template
185 <
186  template<class> class ReactionType,
187  class ReactionThermo,
188  class ReactionRate
189 >
190 Foam::scalar
192 <
193  ReactionType,
194  ReactionThermo,
195  ReactionRate
196 >::kr
197 (
198  const scalar p,
199  const scalar T,
200  const scalarField& c
201 ) const
202 {
203  return rk_(p, T, c);
204 }
205 
206 
207 template
208 <
209  template<class> class ReactionType,
210  class ReactionThermo,
211  class ReactionRate
212 >
214 <
215  ReactionType,
216  ReactionThermo,
217  ReactionRate
218 >::write
219 (
220  Ostream& os
221 ) const
222 {
224 
225  os << indent << "forward" << nl;
226  os << indent << token::BEGIN_BLOCK << nl;
227  os << incrIndent;
228  fk_.write(os);
229  os << decrIndent;
230  os << indent << token::END_BLOCK << nl;
231 
232  os << indent << "reverse" << nl;
233  os << indent << token::BEGIN_BLOCK << nl;
234  os << incrIndent;
235  rk_.write(os);
236  os << decrIndent;
237  os << indent << token::END_BLOCK << nl;
238 }
239 
240 
241 // ************************************************************************* //
dictionary dict
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:223
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
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
Simple extension of Reaction to handle reversible reactions using equilibrium thermodynamics.
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:262
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:237
A wordList with hashed indices for faster lookup by name.
Info<< "Creating reaction model\n"<< endl;autoPtr< combustionModels::psiCombustionModel > reaction(combustionModels::psiCombustionModel::New(mesh))
virtual Ostream & write(const token &)=0
Write next token to stream.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:230
runTime write()