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-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 
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  const dictionary& dict
60 )
61 :
62  ReactionType<ReactionThermo>(species, thermoDatabase, dict),
63  k_(species, dict)
64 {}
65 
66 
67 template
68 <
69  template<class> class ReactionType,
70  class ReactionThermo,
71  class ReactionRate
72 >
75 (
77  const speciesTable& species
78 )
79 :
80  ReactionType<ReactionThermo>(rr, species),
81  k_(rr.k_)
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 template
88 <
89  template<class> class ReactionType,
90  class ReactionThermo,
91  class ReactionRate
92 >
93 Foam::scalar Foam::ReversibleReaction
94 <
95  ReactionType,
96  ReactionThermo,
97  ReactionRate
98 >::kf
99 (
100  const scalar p,
101  const scalar T,
102  const scalarField& c
103 ) const
104 {
105  return k_(p, T, c);
106 }
107 
108 
109 template
110 <
111  template<class> class ReactionType,
112  class ReactionThermo,
113  class ReactionRate
114 >
115 Foam::scalar Foam::ReversibleReaction
116 <
117  ReactionType,
118  ReactionThermo,
119  ReactionRate
120 >::kr
121 (
122  const scalar kfwd,
123  const scalar p,
124  const scalar T,
125  const scalarField& c
126 ) const
127 {
128  scalar Kc = this->Kc(p, T);
129 
130  if (mag(Kc) > VSMALL)
131  {
132  return kfwd/Kc;
133  }
134  else
135  {
136  return 0;
137  }
138 }
139 
140 
141 template
142 <
143  template<class> class ReactionType,
144  class ReactionThermo,
145  class ReactionRate
146 >
147 Foam::scalar Foam::ReversibleReaction
148 <
149  ReactionType,
150  ReactionThermo,
151  ReactionRate
152 >::kr
153 (
154  const scalar p,
155  const scalar T,
156  const scalarField& c
157 ) const
158 {
159  return kr(kf(p, T, c), p, T, c);
160 }
161 
162 
163 template
164 <
165  template<class> class ReactionType,
166  class ReactionThermo,
167  class ReactionRate
168 >
170 <
171  ReactionType,
172  ReactionThermo,
173  ReactionRate
174 >::write
175 (
176  Ostream& os
177 ) const
178 {
180  k_.write(os);
181 }
182 
183 
184 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
ReversibleReaction(const ReactionType< ReactionThermo > &reaction, const ReactionRate &k)
Construct from components.
A HashTable specialization for hashing pointers.
Definition: HashPtrTable.H:50
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
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\"<< endl;autoPtr< combustionModels::psiCombustionModel > reaction(combustionModels::psiCombustionModel::New(mesh))
dimensioned< scalar > mag(const dimensioned< Type > &)
Simple extension of Reaction to handle reversible reactions using equilibrium thermodynamics.