chempy.kinetics package¶
This package collects functions useful for studying chemical kinetics problems.
Submodules¶
chempy.kinetics.analysis module¶
-
chempy.kinetics.analysis.
dominant_reactions_graph
(concs, rate_exprs_cb, rsys, substance_key, linthreshy=1e-09, fname='dominant_reactions_graph.png', relative=False, combine_equilibria=False, **kwargs)[source]¶
-
chempy.kinetics.analysis.
plot_reaction_contributions
(xyp, rsys, rate_exprs_cb, substance_keys=None, varied=None, axes=None, total=False, linthreshy=1e-09, relative=False, xscale='log', yscale='symlog', xlabel='Time', ylabel=None, combine_equilibria=False, selection=slice(None, None, None), unit_registry=None)[source]¶ Plots per reaction contributions to concentration evolution of a substance.
- Parameters
- xyp
pyodesys.results.Result
instance or length 3 tuple or xout,yout,params - resultpyodesys.results.Result
- substance_keystr
- xyp
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__
(self, T[, constants, units, backend])Evaluates the arrhenius equation for a specified state
count
(self, value, /)Return number of occurrences of value.
from_rateconst_at_T
(Ea, T_k[, backend, …])Constructs an instance from a known rate constant at a given temperature.
index
(self, value[, start, stop])Return first index of value.
Ea_over_R
as_RateExpr
equation_as_string
format
from_fit_of_data
html
unicode
-
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__
(self, state[, constants, units, …])See
chempy.arrhenius.arrhenius_equation()
.count
(self, value, /)Return number of occurrences of value.
from_rateconst_at_T
(\*args, \*\*kwargs)Constructs an instance from a known rate constant at a given temperature.
index
(self, value[, start, stop])Return first index of value.
Ea_over_R
as_RateExpr
equation_as_string
format
from_fit_of_data
html
unicode
-
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
- dHfloat
Enthalpy of activation.
- dSfloat
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__
(self, T[, constants, units, backend])Evaluates the Eyring equation for a specified state
count
(self, value, /)Return number of occurrences of value.
index
(self, value[, start, stop])Return first index of value.
as_RateExpr
dH_over_R
equation_as_string
format
kB_h_times_exp_dS_R
-
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__
(self, state[, constants, units, …])See
chempy.eyring.eyring_equation()
.count
(self, value, /)Return number of occurrences of value.
index
(self, value[, start, stop])Return first index of value.
as_RateExpr
dH_over_R
equation_as_string
format
kB_h_times_exp_dS_R
-
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 a few analytic solutions of initial value problems (IVPs) in chemical kinetics (i.e. integrated rate expressions in closed form). The expressions are useful in e.g. regression or for comparison with numerical solution of the corresponding ODE system.
-
chempy.kinetics.integrated.
binary_irrev
(t, kf, prod, major, minor, backend=None)[source]¶ Analytic product transient of a irreversible 2-to-1 reaction.
Product concentration vs time from second order irreversible kinetics.
- Parameters
- tfloat, Symbol or array_like
- kfnumber or Symbol
Forward (bimolecular) rate constant.
- prodnumber or Symbol
Initial concentration of the complex.
- majornumber or Symbol
Initial concentration of the more abundant reactant.
- minornumber or Symbol
Initial concentration of the less abundant reactant.
- backendmodule or str
Default is ‘numpy’, can also be e.g.
sympy
.
-
chempy.kinetics.integrated.
binary_irrev_cstr
(t, k, r, p, fr, fp, fv, n=1, backend=None)[source]¶ Analytic solution for
2 A -> n B
in a CSTR.- Parameters
- tarray_like
- kfloat_like
Rate constant
- rfloat_like
Initial concentration of reactant.
- pfloat_like
Initial concentration of product.
- frfloat_like
Concentration of reactant in feed.
- fpfloat_like
Concentration of product in feed.
- fvfloat_like
Feed rate / tank volume ratio.
- nint
- backendmodule or str
Default is ‘numpy’, can also be e.g.
sympy
.
- Returns
- length-2 tuple
concentrations of reactant and product
-
chempy.kinetics.integrated.
binary_rev
(t, kf, kb, prod, major, minor, backend=None)[source]¶ Analytic product transient of a reversible 2-to-1 reaction.
Product concentration vs time from second order reversible kinetics.
- Parameters
- tfloat, Symbol or array_like
Time.
- kfnumber or Symbol
Forward (bimolecular) rate constant.
- kbnumber or Symbol
Backward (unimolecular) rate constant.
- prodnumber or Symbol
Initial concentration of the complex.
- majornumber or Symbol
Initial concentration of the more abundant reactant.
- minornumber or Symbol
Initial concentration of the less abundant reactant.
- backendmodule or str
Default is ‘numpy’, can also be e.g.
sympy
.
-
chempy.kinetics.integrated.
pseudo_irrev
(t, kf, prod, major, minor, backend=None)[source]¶ Analytic product transient of a irreversible pseudo first order reaction.
Product concentration vs time from pseudo-first order irreversible kinetics.
- Parameters
- tfloat, Symbol or array_like
Time.
- kfnumber or Symbol
Forward (bimolecular) rate constant.
- kbnumber or Symbol
Backward (unimolecular) rate constant.
- prodnumber or Symbol
Initial concentration of the complex.
- majornumber or Symbol
Initial concentration of the more abundant reactant.
- minornumber or Symbol
Initial concentration of the less abundant reactant.
- backendmodule or str
Default is ‘numpy’, can also be e.g.
sympy
.
-
chempy.kinetics.integrated.
pseudo_rev
(t, kf, kb, prod, major, minor, backend=None)[source]¶ Analytic product transient of a reversible pseudo first order reaction.
Product concentration vs time from pseudo-first order reversible kinetics.
- Parameters
- tfloat, Symbol or array_like
Time.
- kfnumber or Symbol
Forward (bimolecular) rate constant.
- kbnumber or Symbol
Backward (unimolecular) rate constant.
- prodnumber or Symbol
Initial concentration of the complex.
- majornumber or Symbol
Initial concentration of the more abundant reactant.
- minornumber or Symbol
Initial concentration of the less abundant reactant.
- backendmodule or str
Default is ‘numpy’, can also be e.g.
sympy
.
-
chempy.kinetics.integrated.
unary_irrev_cstr
(t, k, r, p, fr, fp, fv, backend=None)[source]¶ Analytic solution for
A -> B
in a CSTR.Analytic solution for a first order process in a continuously stirred tank reactor (CSTR).
- Parameters
- tarray_like
- kfloat_like
Rate constant
- rfloat_like
Initial concentration of reactant.
- pfloat_like
Initial concentration of product.
- frfloat_like
Concentration of reactant in feed.
- fpfloat_like
Concentration of product in feed.
- fvfloat_like
Feed rate / tank volume ratio.
- backendmodule or str
Default is ‘numpy’, can also be e.g.
sympy
.
- Returns
- length-2 tuple
concentrations of reactant and product
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.
chained_parameter_variation
(*args, **kwargs)[source]¶ chained_parameter_variation is deprecated since (not including) 0.5.3, it will be missing in 0.8.0. Use pyodesys.chained_parameter_variation instead.
Integrate an ODE-system for a serie of durations with some parameters changed in-between
- Parameters
- odesys
pyodesys.ODESys
instance durations : iterable of floats init_conc : dict or array_like varied_params : dict mapping parameter name to array_like
Each array_like need to be of same length as durations.
- default_paramsdict or array_like
Default values for the parameters of the ODE system.
- integrate_kwargsdict
Keyword arguments passed on to
pyodesys.ODESys.integrate()
.
- odesys
-
chempy.kinetics.ode.
dCdt_list
(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_list(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, cstr=False, constants=None, **kwargs)[source]¶ Creates a
pyneqsys.SymbolicSys
from aReactionSystem
The parameters passed to RateExpr will contain the key
'time'
corresponding to the independent variable of the IVP.- Parameters
- rsysReactionSystem
Each reaction of
rsys
will have theirReaction.rate_expr()
invoked. Note that ifReaction.param
is not aRateExpr
(or convertible to one throughas_RateExpr()
) it will be used to construct aMassAction
instance.- include_paramsbool (default: True)
Whether rate constants should be included into the rate expressions or left as free parameters in the
pyneqsys.SymbolicSys
instance.- substitutionsdict, optional
Variable substitutions used by rate expressions (in respective Reaction.param). values are allowed to be values of instances of
Expr
.- SymbolicSysclass (optional)
Default :
pyneqsys.SymbolicSys
.- unit_registry: dict (optional)
See
chempy.units.get_derived_units()
.- output_conc_unitunit (Optional)
- output_time_unitunit (Optional)
- cstrbool
Generate expressions for continuously stirred tank reactor.
- constantsmodule
e.g.
chempy.units.default_constants
, parameter keys not found in substitutions will be looked for as an attribute ofconstants
when provided.- **kwargs :
Keyword arguemnts passed on to SymbolicSys.
- Returns
- pyodesys.symbolic.SymbolicSys
- extradict, with keys:
param_keys : list of str instances
unique : OrderedDict mapping str to value (possibly None)
p_units : list of units
max_euler_step_cb : callable or None
linear_dependencies : None or factory of solver callback
rate_exprs_cb : callable
cstr_fr_fc : None or (feed-ratio-key, subtance-key-to-feed-conc-key-map)
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, extra = get_odesys(rsys) >>> 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
- concarray_like
concentrations (floats or symbolic objects)
- rsysReactionSystem instance
See
ReactionSystem
- variablesdict (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 Arrhenius, MassAction >>> rxn.param = MassAction(Arrhenius({'A': 1e10, 'Ea_over_R': 9314})) >>> 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.
Arrhenius
(args=None, unique_keys=None)[source]¶ Bases:
chempy.util._expr.Expr
Rate expression for a Arrhenius-type of rate: c0*exp(-c1/T)
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 = MassAction(Arrhenius([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
- argument_defaults
- nargs
- trivially_zero
Methods
__call__
(self, variables[, backend])Call self as a function.
arg
(self, variables, index[, backend, evaluate])- Parameters
args_dimensionality
(self, reaction)return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
dedimensionalisation
(self, unit_registry[, …])Create an instance with consistent units from a unit_registry
fk
(\*args)Alternative constructor “from keys”, *args is used as
unique_keys
.from_callback
(callback[, attr])Factory of subclasses
latex
(self[, variables, backend, default])- Parameters
all_args
all_parameter_keys
all_params
all_unique_keys
string
-
args_dimensionality
(self, reaction)[source]¶ return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
-
argument_names
= ('A', 'Ea_over_R')¶
-
parameter_keys
= ('temperature',)¶
-
class
chempy.kinetics.rates.
Eyring
(args=None, unique_keys=None)[source]¶ Bases:
chempy.util._expr.Expr
Rate expression for Eyring eq: c0*T*exp(-c1/T)
Note that choice of standard state (c^0) will matter for order > 1.
- Attributes
- nargs
- trivially_zero
Methods
__call__
(self, variables[, backend])Call self as a function.
arg
(self, variables, index[, backend, evaluate])- Parameters
args_dimensionality
(self, reaction)return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
dedimensionalisation
(self, unit_registry[, …])Create an instance with consistent units from a unit_registry
fk
(\*args)Alternative constructor “from keys”, *args is used as
unique_keys
.from_callback
(callback[, attr])Factory of subclasses
latex
(self[, variables, backend, default])- Parameters
all_args
all_parameter_keys
all_params
all_unique_keys
string
-
args_dimensionality
(self, reaction)[source]¶ return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
-
argument_defaults
= (array(1.) * M,)¶
-
argument_names
= ('kB_h_times_exp_dS_R', 'dH_over_R', 'conc0')¶
-
parameter_keys
= ('temperature',)¶
-
class
chempy.kinetics.rates.
EyringHS
(args=None, unique_keys=None)[source]¶ Bases:
chempy.util._expr.Expr
- Attributes
- nargs
- trivially_zero
Methods
__call__
(self, variables[, backend, reaction])Call self as a function.
arg
(self, variables, index[, backend, evaluate])- Parameters
args_dimensionality
(self, \*\*kwargs)return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
dedimensionalisation
(self, unit_registry[, …])Create an instance with consistent units from a unit_registry
fk
(\*args)Alternative constructor “from keys”, *args is used as
unique_keys
.from_callback
(callback[, attr])Factory of subclasses
latex
(self[, variables, backend, default])- Parameters
all_args
all_parameter_keys
all_params
all_unique_keys
string
-
args_dimensionality
(self, **kwargs)[source]¶ return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
-
argument_defaults
= (array(1.) * M,)¶
-
argument_names
= ('dH', 'dS', 'c0')¶
-
parameter_keys
= ('temperature', 'molar_gas_constant', 'Boltzmann_constant', 'Planck_constant')¶
-
class
chempy.kinetics.rates.
MassAction
(args=None, unique_keys=None)[source]¶ Bases:
chempy.kinetics.rates.RateExpr
,chempy.util._expr.UnaryWrapper
Rate-expression of mass-action type
Notes
__call__()
requires aReaction
instance to be passed asreaction
keyword argument.Examples
>>> ma = MassAction([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
- argument_defaults
- nargs
- trivially_zero
Methods
__call__
(self, variables[, backend, reaction])Call self as a function.
arg
(self, variables, index[, backend, evaluate])- Parameters
args_dimensionality
(self, reaction)return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
dedimensionalisation
(self, unit_registry[, …])Create an instance with consistent units from a unit_registry
fk
(\*args)Alternative constructor “from keys”, *args is used as
unique_keys
.from_callback
(callback[, attr])Factory of subclasses
latex
(self[, variables, backend, default])- Parameters
subclass_from_callback
(\*args, \*\*kwargs)subclass_from_callback is deprecated.
active_conc_prod
all_args
all_parameter_keys
all_params
all_unique_keys
get_named_keys
rate_coeff
string
-
active_conc_prod
(self, variables, backend=<module 'math' from '/opt/cpython-3.7/lib/python3.7/lib-dynload/math.cpython-37m-x86_64-linux-gnu.so'>, reaction=None)[source]¶
-
args_dimensionality
(self, reaction)[source]¶ return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
-
argument_names
= ('rate_constant',)¶
-
classmethod
from_callback
(callback, attr='__call__', **kwargs)[source]¶ Factory of subclasses
- Parameters
- callbackcallable
signature: *args, backend=math
- attrstr
What attribute to override
- argument_namestuple of str, optional
- argument_defaultstuple of floats, optional
- parameter_keystuple of str, optional,
- nargsint, optional
Examples
>>> from operator import add; from functools import reduce >>> def poly(args, x, backend=math): ... x0 = args[0] ... return reduce(add, [c*(x-x0)**i for i, c in enumerate(args[1:])]) ... >>> Poly = Expr.from_callback(poly, parameter_keys=('x',), argument_names=('x0', Ellipsis)) >>> p = Poly([1, 3, 2, 5]) >>> p({'x': 7}) == 3 + 2*(7-1) + 5*(7-1)**2 True >>> q = Poly([1, 3, 2, 5], unique_keys=('x0_q',)) >>> q({'x': 7, 'x0_q': 0}) == 3 + 2*7 + 5*7**2 True
-
class
chempy.kinetics.rates.
Radiolytic
(args=None, unique_keys=None)¶ Bases:
chempy.kinetics.rates.RadiolyticBase
- Attributes
- argument_defaults
- nargs
- trivially_zero
Methods
arg
(self, variables, index[, backend, evaluate])- Parameters
dedimensionalisation
(self, unit_registry[, …])Create an instance with consistent units from a unit_registry
fk
(\*args)Alternative constructor “from keys”, *args is used as
unique_keys
.from_callback
(callback[, attr])Factory of subclasses
g_value
(\*args, \*\*kwargs)g_value is deprecated.
latex
(self[, variables, backend, default])- Parameters
subclass_from_callback
(\*args, \*\*kwargs)subclass_from_callback is deprecated.
__call__
all_args
all_parameter_keys
all_params
all_unique_keys
args_dimensionality
g_values
string
-
args_dimensionality
(self, reaction)¶
-
argument_names
= ('radiolytic_yield',)¶
-
g_value
(*args, **kwargs)¶ g_value is deprecated. Use Radiolytic.all_args instead.
-
g_values
(self, *args, **kwargs)¶
-
parameter_keys
= ('density', 'doserate')¶
-
class
chempy.kinetics.rates.
RadiolyticBase
(args=None, unique_keys=None)[source]¶ Bases:
chempy.kinetics.rates.RateExpr
- Attributes
- argument_defaults
- argument_names
- nargs
- trivially_zero
Methods
__call__
(self, variables[, backend])Call self as a function.
arg
(self, variables, index[, backend, evaluate])- Parameters
args_dimensionality
(self, \*\*kwargs)return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
dedimensionalisation
(self, unit_registry[, …])Create an instance with consistent units from a unit_registry
fk
(\*args)Alternative constructor “from keys”, *args is used as
unique_keys
.from_callback
(callback[, attr])Factory of subclasses
latex
(self[, variables, backend, default])- Parameters
subclass_from_callback
(\*args, \*\*kwargs)subclass_from_callback is deprecated.
all_args
all_parameter_keys
all_params
all_unique_keys
string
-
class
chempy.kinetics.rates.
RampedTemp
(args=None, unique_keys=None)[source]¶ Bases:
chempy.util._expr.Expr
Ramped temperature, pass as substitution to e.g.
get_odesys
- Attributes
- argument_defaults
- nargs
- trivially_zero
Methods
__call__
(self, variables[, backend])Call self as a function.
arg
(self, variables, index[, backend, evaluate])- Parameters
args_dimensionality
(self, \*\*kwargs)return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
dedimensionalisation
(self, unit_registry[, …])Create an instance with consistent units from a unit_registry
fk
(\*args)Alternative constructor “from keys”, *args is used as
unique_keys
.from_callback
(callback[, attr])Factory of subclasses
latex
(self[, variables, backend, default])- Parameters
all_args
all_parameter_keys
all_params
all_unique_keys
string
-
args_dimensionality
(self, **kwargs)[source]¶ return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
-
argument_names
= ('T0', 'dTdt')¶
-
parameter_keys
= ('time',)¶
-
class
chempy.kinetics.rates.
RateExpr
(args=None, unique_keys=None)[source]¶ Bases:
chempy.util._expr.Expr
Baseclass for rate expressions, see source code of e.g. MassAction & Radiolytic.
- Attributes
- argument_defaults
- argument_names
- nargs
- trivially_zero
Methods
__call__
(self, variables[, backend])Call self as a function.
arg
(self, variables, index[, backend, evaluate])- Parameters
args_dimensionality
(self, \*\*kwargs)return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
dedimensionalisation
(self, unit_registry[, …])Create an instance with consistent units from a unit_registry
fk
(\*args)Alternative constructor “from keys”, *args is used as
unique_keys
.from_callback
(callback[, attr])Factory of subclasses
latex
(self[, variables, backend, default])- Parameters
subclass_from_callback
(\*args, \*\*kwargs)subclass_from_callback is deprecated.
all_args
all_parameter_keys
all_params
all_unique_keys
string
-
classmethod
subclass_from_callback
(*args, **kwargs)[source]¶ subclass_from_callback is deprecated. Use from_callback instead.
Override RateExpr.__call__
- Parameters
- cbcallback
With signature (variables, all_args, backend) -> scalar where variables is a dict, all_args a tuple and backend a module.
- cls_attrsdict, 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]) >>> print('%.5g' % k({'temperature': 298.15, 'O2': 1.1e-3, 'rxn': rxn})) 22.186
-
class
chempy.kinetics.rates.
SinTemp
(args=None, unique_keys=None)[source]¶ Bases:
chempy.util._expr.Expr
- Attributes
- argument_defaults
- nargs
- trivially_zero
Methods
__call__
(self, variables[, backend])Call self as a function.
arg
(self, variables, index[, backend, evaluate])- Parameters
args_dimensionality
(self, \*\*kwargs)return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
dedimensionalisation
(self, unit_registry[, …])Create an instance with consistent units from a unit_registry
fk
(\*args)Alternative constructor “from keys”, *args is used as
unique_keys
.from_callback
(callback[, attr])Factory of subclasses
latex
(self[, variables, backend, default])- Parameters
all_args
all_parameter_keys
all_params
all_unique_keys
string
-
args_dimensionality
(self, **kwargs)[source]¶ return tuple of dicts mapping str to int (‘length’, ‘mass’, ‘time’, ‘current’, ‘temperature’, ‘luminous_intensity’, ‘amount’)
-
argument_names
= ('Tbase', 'Tamp', 'angvel', 'phase')¶
-
parameter_keys
= ('time',)¶
-
chempy.kinetics.rates.
mk_Radiolytic
(*doserate_names)[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.
- Parameters
- *doserate_namesstr instances
Default: (‘’,)
Notes
The instance __call__ will require by default
'density'
and'doserate'
in variables.Examples
>>> RadiolyticAlpha = mk_Radiolytic('alpha') >>> RadiolyticGamma = mk_Radiolytic('gamma') >>> dihydrogen_alpha = RadiolyticAlpha([0.8e-7]) >>> dihydrogen_gamma = RadiolyticGamma([0.45e-7]) >>> RadiolyticAB = mk_Radiolytic('alpha', 'beta')