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-2024 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 ThermoType, class ReactionRate>
33 (
35  const ReactionRate& forwardReactionRate,
36  const ReactionRate& reverseReactionRate
37 )
38 :
39  Reaction<ThermoType>(reaction),
40  kf_(forwardReactionRate),
41  kr_(reverseReactionRate)
42 {}
43 
44 
45 template<class ThermoType, class ReactionRate>
48 (
49  const speciesTable& species,
50  const PtrList<ThermoType>& speciesThermo,
51  const dictionary& dict
52 )
53 :
54  Reaction<ThermoType>(species, speciesThermo, dict),
55  kf_(species, this->kfDims(), dict.subDict("forward")),
56  kr_(species, this->krDims(), dict.subDict("reverse"))
57 {}
58 
59 
60 template<class ThermoType, class ReactionRate>
63 (
64  const speciesTable& species,
65  const PtrList<ThermoType>& speciesThermo,
66  const objectRegistry& ob,
67  const dictionary& dict
68 )
69 :
70  Reaction<ThermoType>(species, speciesThermo, dict),
71  kf_(species, ob, this->kfDims(), dict.subDict("forward")),
72  kr_(species, ob, this->krDims(), dict.subDict("reverse"))
73 {}
74 
75 
76 template<class ThermoType, class ReactionRate>
79 (
81  nerr,
82  const speciesTable& species
83 )
84 :
85  Reaction<ThermoType>(nerr, species),
86  kf_(nerr.kf_),
87  kr_(nerr.kr_)
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
93 template<class ThermoType, class ReactionRate>
95 preEvaluate() const
96 {
97  kf_.preEvaluate();
98  kr_.preEvaluate();
99 }
100 
101 
102 template<class ThermoType, class ReactionRate>
104 postEvaluate() const
105 {
106  kf_.postEvaluate();
107  kr_.postEvaluate();
108 }
109 
110 
111 template<class ThermoType, class ReactionRate>
112 Foam::scalar
114 (
115  const scalar p,
116  const scalar T,
117  const scalarField& c,
118  const label li
119 ) const
120 {
121  return kf_(p, T, c, li);
122 }
123 
124 
125 template<class ThermoType, class ReactionRate>
126 Foam::scalar
128 (
129  const scalar,
130  const scalar p,
131  const scalar T,
132  const scalarField& c,
133  const label li
134 ) const
135 {
136  return kr_(p, T, c, li);
137 }
138 
139 
140 template<class ThermoType, class ReactionRate>
141 Foam::scalar
143 (
144  const scalar p,
145  const scalar T,
146  const scalarField& c,
147  const label li
148 ) const
149 {
150  return kr_(p, T, c, li);
151 }
152 
153 
154 template<class ThermoType, class ReactionRate>
155 Foam::scalar
157 (
158  const scalar p,
159  const scalar T,
160  const scalarField& c,
161  const label li
162 ) const
163 {
164  return kf_.ddT(p, T, c, li);
165 }
166 
167 
168 template<class ThermoType, class ReactionRate>
169 Foam::scalar
171 (
172  const scalar p,
173  const scalar T,
174  const scalarField& c,
175  const label li,
176  const scalar dkfdT,
177  const scalar kr
178 ) const
179 {
180  return kr_.ddT(p, T, c, li);
181 }
182 
183 
184 template<class ThermoType, class ReactionRate>
186 hasDkdc() const
187 {
188  return kf_.hasDdc() || kr_.hasDdc();
189 }
190 
191 
192 template<class ThermoType, class ReactionRate>
194 (
195  const scalar p,
196  const scalar T,
197  const scalarField& c,
198  const label li,
199  scalarField& dkfdc
200 ) const
201 {
202  kf_.ddc(p, T, c, li, dkfdc);
203 }
204 
205 
206 template<class ThermoType, class ReactionRate>
208 (
209  const scalar p,
210  const scalar T,
211  const scalarField& c,
212  const label li,
213  const scalarField& dkfdc,
214  const scalar kr,
215  scalarField& dkrdc
216 ) const
217 {
218  kr_.ddc(p, T, c, li, dkrdc);
219 }
220 
221 
222 template<class ThermoType, class ReactionRate>
224 (
225  Ostream& os
226 ) const
227 {
229 
230  os << indent << "forward" << nl;
231  os << indent << token::BEGIN_BLOCK << nl;
232  os << incrIndent;
233  kf_.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  kr_.write(os);
241  os << decrIndent;
242  os << indent << token::END_BLOCK << nl;
243 }
244 
245 
246 // ************************************************************************* //
Extension of Reaction to handle non-equilibrium reversible reactions.
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 forward rate constant.
virtual scalar dkfdT(const scalar p, const scalar T, const scalarField &c, const label li) const
Temperature derivative of forward rate.
virtual void preEvaluate() const
Pre-evaluation hook.
virtual void postEvaluate() const
Post-evaluation hook.
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.
void dkfdc(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dkfdc) const
Concentration derivative of forward rate.
virtual bool hasDkdc() const
Does this reaction have concentration-dependent rate constants?
NonEquilibriumReversibleReaction(const Reaction< ThermoType > &reaction, const ReactionRate &forwardReactionRate, const ReactionRate &reverseReactionRate)
Construct from components.
virtual scalar kf(const scalar p, const scalar T, const scalarField &c, const label li) const
Forward rate constant.
void dkrdc(const scalar p, const scalar T, const scalarField &c, const label li, const scalarField &dkfdc, const scalar kr, scalarField &dkrdc) const
Concentration derivative of reverse rate.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual Ostream & write(const char)=0
Write character.
Simple extension of ThermoType to handle reaction kinetics in addition to the equilibrium thermodynam...
Definition: Reaction.H:72
virtual void write(Ostream &) const
Write.
Definition: Reaction.C:550
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
A wordList with hashed indices for faster lookup by name.
Registry of regIOobjects.
Reaction base-class holding the specie names and coefficients.
Definition: reaction.H:57
@ BEGIN_BLOCK
Definition: token.H:113
@ END_BLOCK
Definition: token.H:114
const dimensionedScalar c
Speed of light in a vacuum.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:241
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 & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:234
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p