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-2017 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 template
56 <
57  template<class> class ReactionType,
58  class ReactionThermo,
59  class ReactionRate
60 >
62 <
63  ReactionType,
64  ReactionThermo,
65  ReactionRate
66 >::
68 (
69  const speciesTable& species,
70  const HashPtrTable<ReactionThermo>& thermoDatabase,
71  const dictionary& dict
72 )
73 :
74  ReactionType<ReactionThermo>(species, thermoDatabase, dict),
75  fk_(species, dict.subDict("forward")),
76  rk_(species, dict.subDict("reverse"))
77 {}
78 
79 
80 template
81 <
82  template<class> class ReactionType,
83  class ReactionThermo,
84  class ReactionRate
85 >
87 <
88  ReactionType,
89  ReactionThermo,
90  ReactionRate
91 >::
93 (
95  <
96  ReactionType,
97  ReactionThermo,
98  ReactionRate
99  >& nerr,
100  const speciesTable& species
101 )
102 :
103  ReactionType<ReactionThermo>(nerr, species),
104  fk_(nerr.fk_),
105  rk_(nerr.rk_)
106 {}
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
111 template
112 <
113  template<class> class ReactionType,
114  class ReactionThermo,
115  class ReactionRate
116 >
117 Foam::scalar
119 <
120  ReactionType,
121  ReactionThermo,
122  ReactionRate
123 >::kf
124 (
125  const scalar p,
126  const scalar T,
127  const scalarField& c
128 ) const
129 {
130  return fk_(p, T, c);
131 }
132 
133 
134 template
135 <
136  template<class> class ReactionType,
137  class ReactionThermo,
138  class ReactionRate
139 >
140 Foam::scalar
142 <
143  ReactionType,
144  ReactionThermo,
145  ReactionRate
146 >::kr
147 (
148  const scalar,
149  const scalar p,
150  const scalar T,
151  const scalarField& c
152 ) const
153 {
154  return rk_(p, T, c);
155 }
156 
157 
158 template
159 <
160  template<class> class ReactionType,
161  class ReactionThermo,
162  class ReactionRate
163 >
164 Foam::scalar
166 <
167  ReactionType,
168  ReactionThermo,
169  ReactionRate
170 >::kr
171 (
172  const scalar p,
173  const scalar T,
174  const scalarField& c
175 ) const
176 {
177  return rk_(p, T, c);
178 }
179 
180 
181 template
182 <
183  template<class> class ReactionType,
184  class ReactionThermo,
185  class ReactionRate
186 >
188 <
189  ReactionType,
190  ReactionThermo,
191  ReactionRate
192 >::write
193 (
194  Ostream& os
195 ) const
196 {
198 
199  os << indent << "forward" << nl;
200  os << indent << token::BEGIN_BLOCK << nl;
201  os << incrIndent;
202  fk_.write(os);
203  os << decrIndent;
204  os << indent << token::END_BLOCK << nl;
205 
206  os << indent << "reverse" << nl;
207  os << indent << token::BEGIN_BLOCK << nl;
208  os << incrIndent;
209  rk_.write(os);
210  os << decrIndent;
211  os << indent << token::END_BLOCK << nl;
212 }
213 
214 
215 // ************************************************************************* //
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
A HashTable specialization for hashing pointers.
Definition: HashPtrTable.H:50
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
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\"<< 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