chempy.kinetics package¶
This package collects functions useful for studying chemical kinetics problems.
Submodules¶
chempy.kinetics.arrhenius module¶
Contains functions for the Arrhenius equation
(arrhenius_equation()
) and a convenience fitting routine
(fit_arrhenius_equation()
).
-
class
chempy.kinetics.arrhenius.
ArrheniusParam
[source]¶ Bases:
chempy.util.pyutil.ArrheniusParam
Kinetic data in the form of an Arrhenius parameterisation
Parameters: Ea: float
activation energy
A: float
preexponential prefactor (Arrhenius type eq.)
ref: object (default: None)
arbitrary reference (e.g. string representing citation key)
Examples
>>> k = ArrheniusParam(1e13, 40e3) >>> '%.5g' % k(298.15) '9.8245e+05'
Attributes
A
Alias for field number 0 Ea
Alias for field number 1 ref
Alias for field number 2 Methods
__call__
(T[, constants, units, backend])Evaluates the arrhenius equation for a specified state count
(...)equation_as_string
(precision[, tex])format
(precision[, tex])from_fit_of_data
(T, k[, kerr])index
((value, [start, ...)Raises ValueError if the value is not present.
-
class
chempy.kinetics.arrhenius.
ArrheniusParamWithUnits
[source]¶ Bases:
chempy.kinetics.arrhenius.ArrheniusParam
Attributes
A
Alias for field number 0 Ea
Alias for field number 1 ref
Alias for field number 2 Methods
__call__
(state[, constants, units, backend])See chempy.arrhenius.arrhenius_equation()
.count
(...)equation_as_string
(precision[, tex])format
(precision[, tex])from_fit_of_data
(T, k[, kerr])index
((value, [start, ...)Raises ValueError if the value is not present.
-
chempy.kinetics.arrhenius.
arrhenius_equation
(A, Ea, T, constants=None, units=None, backend=None)[source]¶ Returns the rate coefficient according to the Arrhenius equation
Parameters: A: float with unit
frequency factor
Ea: float with unit
activation energy
T: float with unit
temperature
constants: object (optional, default: None)
- if None:
T assumed to be in Kelvin, Ea in J/(K mol)
- else:
attributes accessed: molar_gas_constant Tip: pass quantities.constants
units: object (optional, default: None)
attributes accessed: Joule, Kelvin and mol
backend: module (optional)
module with “exp”, default: numpy, math
chempy.kinetics.eyring module¶
Contains functions for the Eyring equation.
-
class
chempy.kinetics.eyring.
EyringParam
[source]¶ Bases:
chempy.util.pyutil.EyringParam
Kinetic data in the form of an Eyring parameterisation
Parameters: dH : float
Enthalpy of activation.
dS : float
Entropy of activation.
ref: object (default: None)
arbitrary reference (e.g. citation key or dict with bibtex entries)
Examples
>>> k = EyringParam(72e3, 61.4) >>> '%.5g' % k(298.15) '2435.4'
Attributes
dH
Alias for field number 0 dS
Alias for field number 1 ref
Alias for field number 2 Methods
__call__
(T[, constants, units, backend])Evaluates the Eyring equation for a specified state count
(...)equation_as_string
(precision[, tex])format
(precision[, tex])index
((value, [start, ...)Raises ValueError if the value is not present.
-
class
chempy.kinetics.eyring.
EyringParamWithUnits
[source]¶ Bases:
chempy.kinetics.eyring.EyringParam
Attributes
dH
Alias for field number 0 dS
Alias for field number 1 ref
Alias for field number 2 Methods
__call__
(state[, constants, units, backend])See chempy.eyring.eyring_equation()
.count
(...)equation_as_string
(precision[, tex])format
(precision[, tex])index
((value, [start, ...)Raises ValueError if the value is not present.
-
chempy.kinetics.eyring.
eyring_equation
(dH, dS, T, constants=None, units=None, backend=None)[source]¶ Returns the rate coefficient according to the Eyring equation
Parameters: dH: float with unit
Enthalpy of activation.
dS: float with unit
Entropy of activation.
T: float with unit
temperature
constants: object (optional, default: None)
- if None:
T assumed to be in Kelvin, Ea in J/(K mol)
- else:
attributes accessed: molar_gas_constant Tip: pass quantities.constants
units: object (optional, default: None)
attributes accessed: Joule, Kelvin and mol
backend: module (optional)
module with “exp”, default: numpy, math
chempy.kinetics.integrated module¶
This module collects special cases of systems for which analytic solutions are available.
-
chempy.kinetics.integrated.
binary_irrev
(t, kf, P0, t0, excess_C, limiting_C, eps_l, backend=None)[source]¶
-
chempy.kinetics.integrated.
binary_rev
(t, kf, P0, t0, excess_C, limiting_C, eps_l, beta, backend=None)[source]¶
chempy.kinetics.ode module¶
This module contains functions for formulating systems of Ordinary Differential Equations (ODE-systems) which may be integrated numerically to model temporal evolution of concentrations in reaction systems.
-
chempy.kinetics.ode.
dCdt
(rsys, rates)[source]¶ Returns a list of the time derivatives of the concentrations
Parameters: rsys: ReactionSystem instance
rates: array_like
rates (to be weighted by stoichiometries) of the reactions in
rsys
Examples
>>> from chempy import ReactionSystem, Reaction >>> line, keys = 'H2O -> H+ + OH- ; 1e-4', 'H2O H+ OH-' >>> rsys = ReactionSystem([Reaction.from_string(line, keys)], keys) >>> dCdt(rsys, [0.0054]) [-0.0054, 0.0054, 0.0054]
-
chempy.kinetics.ode.
get_odesys
(rsys, include_params=True, substitutions=None, SymbolicSys=None, unit_registry=None, output_conc_unit=None, output_time_unit=None, **kwargs)[source]¶ Creates a
pyneqsys.SymbolicSys
from aReactionSystem
Parameters: rsys : ReactionSystem
note that if
param
if not RateExpr (or convertible to one through_as_RateExpr()
) it will be used to construct aMassAction
instance.include_params : bool (default: False)
whether rate constants should be included into the rate expressions or left as free parameters in the
pyneqsys.SymbolicSys
instance.substitutions : dict, optional
variable substitutions used by rate expressions (in respective Reaction.param). values are allowed to be tuple like: (new_vars, callback)
SymbolicSys : class (optional)
default :
pyneqsys.SymbolicSys
unit_registry: dict (optional)
see
chempy.units.get_derived_units()
output_conc_unit : unit (Optional)
output_time_unit : unit (Optional)
**kwargs :
Keyword arguemnts pass on to SymbolicSys
Returns: pyodesys.symbolic.SymbolicSys
param_keys
unique_keys
p_units
Examples
>>> from chempy import Equilibrium, ReactionSystem >>> eq = Equilibrium({'Fe+3', 'SCN-'}, {'FeSCN+2'}, 10**2) >>> substances = 'Fe+3 SCN- FeSCN+2'.split() >>> rsys = ReactionSystem(eq.as_reactions(kf=3.0), substances) >>> odesys = get_odesys(rsys)[0] >>> init_conc = {'Fe+3': 1.0, 'SCN-': .3, 'FeSCN+2': 0} >>> tout, Cout, info = odesys.integrate(5, init_conc) >>> Cout[-1, :].round(4) array([ 0.7042, 0.0042, 0.2958])
-
chempy.kinetics.ode.
law_of_mass_action_rates
(conc, rsys, variables=None)[source]¶ Returns a generator of reaction rate expressions
Rates from the law of mass action (
Reaction.inact_reac
ignored) from aReactionSystem
.Parameters: conc : array_like
concentrations (floats or symbolic objects)
rsys : ReactionSystem instance
See
ReactionSystem
variables : dict (optional)
to override parameters in the rate expressions of the reactions
Examples
>>> from chempy import ReactionSystem, Reaction >>> line, keys = 'H2O -> H+ + OH- ; 1e-4', 'H2O H+ OH-' >>> rxn = Reaction.from_string(line, keys) >>> rsys = ReactionSystem([rxn], keys) >>> next(law_of_mass_action_rates([55.4, 1e-7, 1e-7], rsys)) 0.00554 >>> from chempy.kinetics.rates import ArrheniusMassAction >>> rxn.param = ArrheniusMassAction({'A': 1e10, 'Ea_over_R': 9314}, rxn=rxn) >>> print('%.5g' % next(law_of_mass_action_rates([55.4, 1e-7, 1e-7], rsys, {'temperature': 293}))) 0.0086693
chempy.kinetics.rates module¶
This module collects object representing rate expressions. It is based
on the chemp.util._expr
module. The API is somewhat cumbersome since
it tries to be compatible with pure python, SymPy and the underlying
units library of ChemPy (quantities
). Consider the API to be provisional.
-
class
chempy.kinetics.rates.
ArrheniusMassAction
(args, unique_keys=None, **kwargs)[source]¶ Bases:
chempy.kinetics.rates.MassAction
Rate expression for a Arrhenius-type of rate
Examples
>>> from math import exp >>> from chempy import Reaction >>> from chempy.units import allclose, default_units as u >>> A = 1e11 / u.second >>> Ea_over_R = 42e3/8.3145 * u.K**-1 >>> ratex = ArrheniusMassAction([A, Ea_over_R]) >>> rxn = Reaction({'R'}, {'P'}, ratex) >>> dRdt = rxn.rate({'R': 3*u.M, 'temperature': 298.15*u.K})['R'] >>> allclose(dRdt, -3*1e11*exp(-42e3/8.3145/298.15)*u.M/u.s) True
Attributes
kwargs
rxn
nargs print_name Methods
__call__
(variables[, backend])as_mass_action
(variables[, backend])rate_coeff
(variables[, backend])subclass_from_callback
(cb[, cls_attrs])Override MassAction.rate_coeff -
argument_names
= ('A', 'Ea_over_R')¶
-
parameter_keys
= ('temperature',)¶
-
-
class
chempy.kinetics.rates.
EyringMassAction
(args, unique_keys=None, **kwargs)[source]¶ Bases:
chempy.kinetics.rates.ArrheniusMassAction
Attributes
kwargs
rxn
nargs print_name Methods
__call__
(variables[, backend])as_mass_action
(variables[, backend])rate_coeff
(variables[, backend])subclass_from_callback
(cb[, cls_attrs])Override MassAction.rate_coeff -
argument_names
= ('kB_h_times_exp_dS_R', 'dH_over_R')¶
-
-
class
chempy.kinetics.rates.
MassAction
(args, unique_keys=None, **kwargs)[source]¶ Bases:
chempy.kinetics.rates.RateExpr
Rate-expression of mass-action type
Examples
>>> ma = MassAction([3.14]) >>> ma.rate_coeff({}) 3.14 >>> from chempy import Reaction >>> r = Reaction.from_string('3 A -> B', param=ma) >>> r.rate({'A': 2}) == {'A': -75.36, 'B': 25.12} True
Attributes
kwargs
rxn
nargs print_name Methods
__call__
(variables[, backend])as_mass_action
(variables[, backend])rate_coeff
(variables[, backend])subclass_from_callback
(cb[, cls_attrs])Override MassAction.rate_coeff -
argument_names
= ('rate_constant',)¶
-
classmethod
subclass_from_callback
(cb, cls_attrs=None)[source]¶ Override MassAction.rate_coeff
Parameters: cb : callback
With signature (variables, all_args, backend) -> scalar where variables is a dict, all_args a tuple and backend a module.
cls_attrs : dict, optional
Attributes to set in subclass, e.g. parameter_keys, nargs
Examples
>>> from functools import reduce >>> from operator import add >>> from chempy import Reaction # d[H2]/dt = 10**(p0 + p1/T + p2/T**2)*[e-]**2 >>> rxn = Reaction({'e-': 2}, {'OH-': 2, 'H2': 1}, None, {'H2O': 2}) >>> def cb(variables, all_args, backend): ... T = variables['temperature'] ... return 10**reduce(add, [p*T**-i for i, p in enumerate(all_args)]) >>> MyMassAction = MassAction.subclass_from_callback(cb, dict(parameter_keys=('temperature',), nargs=-1)) >>> k = MyMassAction([9, 300, -75000], rxn=rxn) >>> print('%.5g' % k({'temperature': 293., 'e-': 1e-10})) 1.4134e-11
-
-
chempy.kinetics.rates.
Radiolytic
¶ alias of
_Radiolytic
-
class
chempy.kinetics.rates.
RadiolyticBase
(args, unique_keys=None, **kwargs)[source]¶ Bases:
chempy.kinetics.rates.RateExpr
Attributes
kwargs
rxn
argument_names nargs print_name Methods
subclass_from_callback
(*args, **kwargs)subclass_from_callback is deprecated. Use from_callback instead.
-
class
chempy.kinetics.rates.
RateExpr
(args, unique_keys=None, **kwargs)[source]¶ Bases:
chempy.util._expr.Expr
Baseclass for rate expressions, see source code of e.g. MassAction & Radiolytic.
Attributes
kwargs
rxn
argument_names nargs print_name Methods
subclass_from_callback
(*args, **kwargs)subclass_from_callback is deprecated. Use from_callback instead. -
kw
= {'ref': None, 'rxn': None}¶
-
rxn
¶
-
classmethod
subclass_from_callback
(*args, **kwargs)[source]¶ subclass_from_callback is deprecated. Use from_callback instead.
Override RateExpr.__call__
Parameters: cb : callback
With signature (variables, all_args, backend) -> scalar where variables is a dict, all_args a tuple and backend a module.
- cls_attrs : dict, optional
Attributes to set in subclass, e.g. parameter_keys, nargs
Examples
>>> from chempy import Reaction >>> rxn = Reaction({'O2': 1, 'H2': 1}, {'H2O2': 1}) # d[H2O2]/dt = p0*exp(-p1/T)*sqrt([O2]) >>> def cb(variables, all_args, backend): ... O2, T = variables['O2'], variables['temperature'] ... p0, p1 = all_args ... return p0*backend.sqrt(O2)*backend.exp(-p1/T) >>> MyRateExpr = RateExpr.subclass_from_callback(cb, dict(parameter_keys=('temperature',),nargs=2)) >>> k = MyRateExpr([1.3e9, 4317.2], rxn=rxn) >>> print('%.5g' % k({'temperature': 298.15, 'O2': 1.1e-3})) 22.186
-
-
chempy.kinetics.rates.
mk_Radiolytic
(*args)[source]¶ Create a Radiolytic rate expression
Note that there is no mass-action dependence in the resulting class, i.e. the rates does not depend on any concentrations.
Examples
>>> RadiolyticAlpha = mk_Radiolytic('doserate_alpha') >>> RadiolyticGamma = mk_Radiolytic('doserate_gamma') >>> dihydrogen_alpha = RadiolyticAlpha([0.8e-7]) >>> dihydrogen_gamma = RadiolyticGamma([0.45e-7])