NonEquilibriumReversibleReaction.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-2019 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<class ReactionThermo, class ReactionRate>
33 (
34  const Reaction<ReactionThermo>& reaction,
35  const ReactionRate& forwardReactionRate,
36  const ReactionRate& reverseReactionRate
37 )
38 :
40  fk_(forwardReactionRate),
41  rk_(reverseReactionRate)
42 {}
43 
44 
45 template<class ReactionThermo, class ReactionRate>
48 (
49  const speciesTable& species,
50  const HashPtrTable<ReactionThermo>& thermoDatabase,
51  const dictionary& dict
52 )
53 :
54  Reaction<ReactionThermo>(species, thermoDatabase, dict),
55  fk_(species, dict.subDict("forward")),
56  rk_(species, dict.subDict("reverse"))
57 {}
58 
59 
60 template<class ReactionThermo, class ReactionRate>
63 (
64  const speciesTable& species,
65  const HashPtrTable<ReactionThermo>& thermoDatabase,
66  const objectRegistry& ob,
67  const dictionary& dict
68 )
69 :
70  Reaction<ReactionThermo>(species, thermoDatabase, dict),
71  fk_(species, ob, dict.subDict("forward")),
72  rk_(species, ob, dict.subDict("reverse"))
73 {}
74 
75 
76 template<class ReactionThermo, class ReactionRate>
79 (
81  const speciesTable& species
82 )
83 :
84  Reaction<ReactionThermo>(nerr, species),
85  fk_(nerr.fk_),
86  rk_(nerr.rk_)
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
92 template<class ReactionThermo, class ReactionRate>
93 Foam::scalar
95 (
96  const scalar p,
97  const scalar T,
98  const scalarField& c,
99  const label li
100 ) const
101 {
102  return fk_(p, T, c, li);
103 }
104 
105 
106 template<class ReactionThermo, class ReactionRate>
107 Foam::scalar
109 (
110  const scalar,
111  const scalar p,
112  const scalar T,
113  const scalarField& c,
114  const label li
115 ) const
116 {
117  return rk_(p, T, c, li);
118 }
119 
120 
121 template<class ReactionThermo, class ReactionRate>
122 Foam::scalar
124 (
125  const scalar p,
126  const scalar T,
127  const scalarField& c,
128  const label li
129 ) const
130 {
131  return rk_(p, T, c, li);
132 }
133 
134 
135 template<class ReactionThermo, class ReactionRate>
136 Foam::scalar
138 (
139  const scalar p,
140  const scalar T,
141  const scalarField& c,
142  const label li
143 ) const
144 {
145  return fk_.ddT(p, T, c, li);
146 }
147 
148 
149 template<class ReactionThermo, class ReactionRate>
150 Foam::scalar
152 (
153  const scalar p,
154  const scalar T,
155  const scalarField& c,
156  const label li,
157  const scalar dkfdT,
158  const scalar kr
159 ) const
160 {
161  return rk_.ddT(p, T, c, li);
162 }
163 
164 
165 template<class ReactionThermo, class ReactionRate>
168 beta() const
169 {
170  return fk_.beta();
171 }
172 
173 
174 template<class ReactionThermo, class ReactionRate>
175 void
177 (
178  const scalar p,
179  const scalar T,
180  const scalarField& c,
181  const label li,
182  scalarField& dcidc
183 ) const
184 {
185  fk_.dcidc(p, T, c, li, dcidc);
186 }
187 
188 
189 template<class ReactionThermo, class ReactionRate>
190 Foam::scalar
192 (
193  const scalar p,
194  const scalar T,
195  const scalarField& c,
196  const label li
197 ) const
198 {
199  return fk_.dcidT(p, T, c, li);
200 }
201 
202 
203 template<class ReactionThermo, class ReactionRate>
204 void
206 (
207  Ostream& os
208 ) const
209 {
211 
212  os << indent << "forward" << nl;
213  os << indent << token::BEGIN_BLOCK << nl;
214  os << incrIndent;
215  fk_.write(os);
216  os << decrIndent;
217  os << indent << token::END_BLOCK << nl;
218 
219  os << indent << "reverse" << nl;
220  os << indent << token::BEGIN_BLOCK << nl;
221  os << incrIndent;
222  rk_.write(os);
223  os << decrIndent;
224  os << indent << token::END_BLOCK << nl;
225 }
226 
227 
228 // ************************************************************************* //
dictionary dict
virtual Ostream & write(const char)=0
Write character.
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
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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
virtual const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
A HashTable specialization for hashing pointers.
Definition: HashPtrTable.H:50
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:56
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:934
virtual scalar dkfdT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of forward rate.
virtual scalar kf(const scalar p, const scalar T, const scalarField &c, const label li) const
Forward rate constant.
NonEquilibriumReversibleReaction(const Reaction< ReactionThermo > &reaction, const ReactionRate &forwardReactionRate, const ReactionRate &reverseReactionRate)
Construct from components.
Simple extension of Reaction to handle reversible reactions using equilibrium thermodynamics.
CombustionModel< rhoReactionThermo > & reaction
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
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 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.
A wordList with hashed indices for faster lookup by name.
Registry of regIOobjects.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
virtual scalar dcidT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of the pressure dependent term.