30 template<
class ThermoType,
class ReactionRate>
35 const ReactionRate& forwardReactionRate,
36 const ReactionRate& reverseReactionRate
40 kf_(forwardReactionRate),
41 kr_(reverseReactionRate)
45 template<
class ThermoType,
class ReactionRate>
55 kf_(species, this->kfDims(),
dict.subDict(
"forward")),
56 kr_(species, this->krDims(),
dict.subDict(
"reverse"))
60 template<
class ThermoType,
class ReactionRate>
71 kf_(species, ob, this->kfDims(),
dict.subDict(
"forward")),
72 kr_(species, ob, this->krDims(),
dict.subDict(
"reverse"))
76 template<
class ThermoType,
class ReactionRate>
93 template<
class ThermoType,
class ReactionRate>
102 template<
class ThermoType,
class ReactionRate>
111 template<
class ThermoType,
class ReactionRate>
121 return kf_(
p,
T,
c, li);
125 template<
class ThermoType,
class ReactionRate>
136 return kr_(
p,
T,
c, li);
140 template<
class ThermoType,
class ReactionRate>
150 return kr_(
p,
T,
c, li);
154 template<
class ThermoType,
class ReactionRate>
164 return kf_.ddT(
p,
T,
c, li);
168 template<
class ThermoType,
class ReactionRate>
180 return kr_.ddT(
p,
T,
c, li);
184 template<
class ThermoType,
class ReactionRate>
188 return kf_.hasDdc() || kr_.hasDdc();
192 template<
class ThermoType,
class ReactionRate>
202 kf_.ddc(
p,
T,
c, li, dkfdc);
206 template<
class ThermoType,
class ReactionRate>
218 kr_.ddc(
p,
T,
c, li, dkrdc);
222 template<
class ThermoType,
class ReactionRate>
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 write(Ostream &) const
Write.
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,...
virtual Ostream & write(const char)=0
Write character.
Simple extension of ThermoType to handle reaction kinetics in addition to the equilibrium thermodynam...
virtual void write(Ostream &) const
Write.
A list of keyword definitions, which are a keyword followed by any number of values (e....
A wordList with hashed indices for faster lookup by name.
Registry of regIOobjects.
Reaction base-class holding the specie names and coefficients.
const dimensionedScalar c
Speed of light in a vacuum.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Ostream & indent(Ostream &os)
Indent stream.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)