ReversibleReaction.C
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-2021 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<class ReactionThermo, class ReactionRate>
33 (
34  const Reaction<ReactionThermo>& reaction,
35  const ReactionRate& k
36 )
37 :
39  k_(k)
40 {}
41 
42 
43 template<class ReactionThermo, class ReactionRate>
46 (
47  const speciesTable& species,
48  const HashPtrTable<ReactionThermo>& thermoDatabase,
49  const dictionary& dict
50 )
51 :
52  Reaction<ReactionThermo>(species, thermoDatabase, dict),
53  k_(species, dict)
54 {}
55 
56 
57 template<class ReactionThermo, class ReactionRate>
60 (
61  const speciesTable& species,
62  const HashPtrTable<ReactionThermo>& thermoDatabase,
63  const objectRegistry& ob,
64  const dictionary& dict
65 )
66 :
67  Reaction<ReactionThermo>(species, thermoDatabase, dict),
68  k_(species, ob, dict)
69 {}
70 
71 
72 template<class ReactionThermo, class ReactionRate>
75 (
77  const speciesTable& species
78 )
79 :
80  Reaction<ReactionThermo>(rr, species),
81  k_(rr.k_)
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 template<class ReactionThermo, class ReactionRate>
88 void
90 {
91  k_.preEvaluate();
92 }
93 
94 
95 template<class ReactionThermo, class ReactionRate>
96 void
98 {
99  k_.postEvaluate();
100 }
101 
102 
103 template<class ReactionThermo, class ReactionRate>
105 (
106  const scalar p,
107  const scalar T,
108  const scalarField& c,
109  const label li
110 ) const
111 {
112  return k_(p, T, c, li);
113 }
114 
115 
116 template<class ReactionThermo, class ReactionRate>
118 (
119  const scalar kfwd,
120  const scalar p,
121  const scalar T,
122  const scalarField& c,
123  const label li
124 ) const
125 {
126  return kfwd/max(this->Kc(p, T), rootSmall);
127 }
128 
129 
130 template<class ReactionThermo, class ReactionRate>
132 (
133  const scalar p,
134  const scalar T,
135  const scalarField& c,
136  const label li
137 ) const
138 {
139  return kr(kf(p, T, c, li), p, T, c, li);
140 }
141 
142 
143 template<class ReactionThermo, class ReactionRate>
145 (
146  const scalar p,
147  const scalar T,
148  const scalarField& c,
149  const label li
150 ) const
151 {
152  return k_.ddT(p, T, c, li);
153 }
154 
155 
156 template<class ReactionThermo, class ReactionRate>
158 (
159  const scalar p,
160  const scalar T,
161  const scalarField& c,
162  const label li,
163  const scalar dkfdT,
164  const scalar kr
165 ) const
166 {
167  scalar Kc = max(this->Kc(p, T), rootSmall);
168 
169  return dkfdT/Kc - kr*this->dKcdTbyKc(p, T);
170 }
171 
172 
173 template<class ReactionThermo, class ReactionRate>
176 {
177  return k_.beta();
178 }
179 
180 
181 template<class ReactionThermo, class ReactionRate>
183 (
184  const scalar p,
185  const scalar T,
186  const scalarField& c,
187  const label li,
188  scalarField& dcidc
189 ) const
190 {
191  k_.dcidc(p, T, c, li, dcidc);
192 }
193 
194 
195 template<class ReactionThermo, class ReactionRate>
197 (
198  const scalar p,
199  const scalar T,
200  const scalarField& c,
201  const label li
202 ) const
203 {
204  return k_.dcidT(p, T, c, li);
205 }
206 
207 
208 template<class ReactionThermo, class ReactionRate>
210 (
211  Ostream& os
212 ) const
213 {
215  k_.write(os);
216 }
217 
218 
219 // ************************************************************************* //
dictionary dict
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
combustionModel & reaction
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:50
virtual scalar dcidT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of the pressure dependent term.
virtual void dcidc(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:55
virtual void preEvaluate() const
Pre-evaluation hook.
virtual scalar kr(const scalar kfwd, const scalar p, const scalar T, const scalarField &c, const label li) const
Reverse rate constant from the given formard rate constant.
virtual void write(Ostream &) const
Write.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
virtual void postEvaluate() const
Post-evaluation hook.
virtual scalar kf(const scalar p, const scalar T, const scalarField &c, const label li) const
Forward rate constant.
virtual scalar dkfdT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of forward rate.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
A wordList with hashed indices for faster lookup by name.
virtual scalar dkrdT(const scalar p, const scalar T, const scalarField &c, const label li, const scalar dkfdT, const scalar kr) const
Temperature derivative of backward rate.
ReversibleReaction(const Reaction< ReactionThermo > &reaction, const ReactionRate &k)
Construct from components.
Registry of regIOobjects.
Simple extension of Reaction to handle reversible reactions using equilibrium thermodynamics.