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-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 
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>
94 preEvaluate() const
95 {
96  fk_.preEvaluate();
97  rk_.preEvaluate();
98 }
99 
100 
101 template<class ReactionThermo, class ReactionRate>
104 {
105  fk_.postEvaluate();
106  rk_.postEvaluate();
107 }
108 
109 
110 template<class ReactionThermo, class ReactionRate>
111 Foam::scalar
113 (
114  const scalar p,
115  const scalar T,
116  const scalarField& c,
117  const label li
118 ) const
119 {
120  return fk_(p, T, c, li);
121 }
122 
123 
124 template<class ReactionThermo, class ReactionRate>
125 Foam::scalar
127 (
128  const scalar,
129  const scalar p,
130  const scalar T,
131  const scalarField& c,
132  const label li
133 ) const
134 {
135  return rk_(p, T, c, li);
136 }
137 
138 
139 template<class ReactionThermo, class ReactionRate>
140 Foam::scalar
142 (
143  const scalar p,
144  const scalar T,
145  const scalarField& c,
146  const label li
147 ) const
148 {
149  return rk_(p, T, c, li);
150 }
151 
152 
153 template<class ReactionThermo, class ReactionRate>
154 Foam::scalar
156 (
157  const scalar p,
158  const scalar T,
159  const scalarField& c,
160  const label li
161 ) const
162 {
163  return fk_.ddT(p, T, c, li);
164 }
165 
166 
167 template<class ReactionThermo, class ReactionRate>
168 Foam::scalar
170 (
171  const scalar p,
172  const scalar T,
173  const scalarField& c,
174  const label li,
175  const scalar dkfdT,
176  const scalar kr
177 ) const
178 {
179  return rk_.ddT(p, T, c, li);
180 }
181 
182 
183 template<class ReactionThermo, class ReactionRate>
186 beta() const
187 {
188  return fk_.beta();
189 }
190 
191 
192 template<class ReactionThermo, class ReactionRate>
193 void
195 (
196  const scalar p,
197  const scalar T,
198  const scalarField& c,
199  const label li,
200  scalarField& dcidc
201 ) const
202 {
203  fk_.dcidc(p, T, c, li, dcidc);
204 }
205 
206 
207 template<class ReactionThermo, class ReactionRate>
208 Foam::scalar
210 (
211  const scalar p,
212  const scalar T,
213  const scalarField& c,
214  const label li
215 ) const
216 {
217  return fk_.dcidT(p, T, c, li);
218 }
219 
220 
221 template<class ReactionThermo, class ReactionRate>
222 void
224 (
225  Ostream& os
226 ) const
227 {
229 
230  os << indent << "forward" << nl;
231  os << indent << token::BEGIN_BLOCK << nl;
232  os << incrIndent;
233  fk_.write(os);
234  os << decrIndent;
235  os << indent << token::END_BLOCK << nl;
236 
237  os << indent << "reverse" << nl;
238  os << indent << token::BEGIN_BLOCK << nl;
239  os << incrIndent;
240  rk_.write(os);
241  os << decrIndent;
242  os << indent << token::END_BLOCK << nl;
243 }
244 
245 
246 // ************************************************************************* //
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:156
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
virtual void postEvaluate() const
Post-evaluation hook.
virtual const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
A HashTable specialisation 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:55
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:982
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.
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.
virtual void preEvaluate() const
Pre-evaluation hook.
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.