python-enzymegraph

A Python package for generating models of enzymes under the quasi-steady-state assumption (QSSA).

Requirements

Currently tested with only Python 3.4.

Requires sympy.

Installation

pip install git+git://github.com/robjstan/python-enzymegraph.git

References

Gunawardena, J.
A linear framework for time-scale separation in nonlinear biochemical systems.
PLoS ONE (2012)

Gabow, H. N. & Myers, E. W.
Finding all spanning trees of directed and undirected graphs.
SIAM Journal on Computing (1978)

Example

Biochemical example is the Michaelis-Menten function with product inhibition. $$E + S \underset{k_2}{\overset{k_1}{\rightleftharpoons}} ES \overset{k_\text{cat}}{\rightleftharpoons} EP \underset{k_4}{\overset{k_3}{\rightleftharpoons}} E + P$$

Python setup


In [1]:
from enzymegraph import *

In [2]:
from sympy import *
from numpy import linspace
from scipy.integrate import odeint
import matplotlib.pyplot as plt
%matplotlib inline

Symbolic setup

First we need to initialise the varaibles that describe the biochemical species and complexes.

e  = free enzyme
s  = substrate
p  = product
es = enzyme bound to substrate
ep = enzyme bound to product

In [3]:
e, s, p, es, ep = symbols('e s p es ep', positive = True)

Then we need to describe the edges between these species/complexes, labelled with their rates, including the uptake of substrate and product.


In [4]:
k1, k2, k3, k4, kcat = symbols('k_1 k_2 k_3, k_4 k_cat', positive = True)

edges = {(e, es): k1*s, (es, e): k2,
         (e, ep): k3*p, (ep, e): k4,
         (es,ep): kcat, }

Creating the graph object

enzymegraph.enzymegraph(...)

The class constructor enzymegraph.enzymegraph(...) will take as its first argument either:

  • a dictionary with keys=edges and values=edge labels
  • or, a list of edges (edge labels assumed to be equal to 1)

In [5]:
graph = enzymegraph(edges)
graph


Out[5]:
<enzymegraph.enzymegraph.enzymegraph at 0x112c9b208>

In [6]:
edges = list(edges.keys())

enzymegraph(edges)


Out[6]:
<enzymegraph.enzymegraph.enzymegraph at 0x112c5ae10>

enzymegraph.from_matrix(...)

The class constructor enzymegraph.from_matrix(...) will take as its first argument a matrix of edges, with an optional second argument giving vertex names.


In [7]:
m = [[ 0,k1*s, k3*p],
     [k2,   0, kcat],
     [k4,   0,    0]]
v = [e, es, ep]

enzymegraph.from_matrix(m, v).edges


Out[7]:
{(e, es): k_1*s, (ep, e): k_4, (es, e): k_2, (e, ep): k_3*p, (es, ep): k_cat}

ODE outputs

graph.ode_model()

The class method enzymegraph.ode_model() returns the equivalent mass-action (ODE) model of the system.


In [8]:
for var, ode in graph.ode_model().items():
    print("d%s/dt = %s" % (var, ode))


dep/dt = e*k_3*p - ep*k_4 + es*k_cat
de/dt = -e*k_1*s - e*k_3*p + ep*k_4 + es*k_2
des/dt = e*k_1*s - es*k_2 - es*k_cat

graph.ode_function()

The method enzymegraph.ode_function() returns the equivalent mass-action (ODE) model of the system as a python function (suitable for simulation using scipy.odeint(...).

This requires three extra parameters to be given:

  • param_dict, a dictionary giving specific values for the system parameters
  • variables, a list of the model variables in a desired order
  • extra_odes, a dictionary of additional ODEs that complete the system (e.g. equations that describe the production of substrate and product)

In [9]:
param_dict = {k1:1, k2:0.1, k3:2, k4:0.5, kcat:1}
variables = [e, s, p, es, ep]
extra_odes = {s: k2*es - k1*e*s, p: k4*ep - k3*e*p}
ode_func = graph.ode_function(param_dict, variables, extra_odes)
ode_func


Out[9]:
<function enzymegraph.enzymegraph.<lambda>>

In [10]:
concs = zip(*odeint(ode_func, [1,1,0,0,0], linspace(0,10,100)))

for conc, var in zip(concs, variables):
    plt.plot(conc, label=var)
plt.xlabel("time")
plt.ylabel("concentration")
plt.legend()
plt.show()


Spanning trees

graph.spanning_trees()

The method enzymegraph.spanning_trees() enumerates all the possible spanning trees of the graph, returning these as enzymegraph objects.


In [11]:
for span_tree in graph.spanning_trees():
    print(span_tree.edges)


{(es, e): k_2, (e, ep): k_3*p}
{(e, ep): k_3*p, (es, ep): k_cat}
{(e, es): k_1*s, (es, ep): k_cat}
{(ep, e): k_4, (es, e): k_2}
{(ep, e): k_4, (es, ep): k_cat}
{(e, es): k_1*s, (ep, e): k_4}

TikZ output

graph.graph_as_tikz(...)

The method enzymegraph.graph_as_tikz(...) outputs the graph as TikZ commands, suitable for inclusion in a LaTeX document.

This requires extra parameters to be given:

  • vertex_pos, a dictionary giving Tikz position arguments for each of the complexes
  • edge_styles (optional), a dictionary giving a Tikz line argument for any of the edges
  • relabel (optional), whether to relabel the output TikZ nodes (for instance if their names are otherwise complicated (and so may cause TikZ errors)
  • vertex_text_formatter (optional), a lambda that takes the symbolic label of the vertex and returns a formatted string.
  • edge_text_formatter (optional), a lambda that takes the symbolic label of the edge and returns a formatted string.

In [12]:
vertex_pos = {e:  "(0,0)",
              es: "(240:2)",
              ep: "(300:2)",}
edge_styles = {(e, es): "to[out=210,in=90]",
               (e, ep): "to[out=330,in=90]",}

print(graph.graph_as_tikz(vertex_pos, edge_styles, relabel=False))


\tikzset{edge label/.style={font=\tiny, fill=white, inner sep=1}}

\begin{tikzpicture}
	\node (e) at (0,0) {$e$};
	\node (ep) at (300:2) {$ep$};
	\node (es) at (240:2) {$es$};
	\begin{scope}[every node/.style=edge label]
		\draw[->] (e) to[out=210,in=90] node {$k_{1} s$} (es);
		\draw[->] (e) to[out=330,in=90] node {$k_{3} p$} (ep);
		\draw[->] (ep) -- node {$k_{4}$} (e);
		\draw[->] (es) -- node {$k_{2}$} (e);
		\draw[->] (es) -- node {$k_{cat}$} (ep);
	\end{scope}
\end{tikzpicture}

graph.spanning_trees_as_tikz(...) (beta)

The method enzymegraph.spanning_trees_as_tikz(...) outputs all the spanning trees of the system as TikZ commands, suitable for inclusion in a LaTeX document.

This requires extra parameters, as described under enzymegraph.graph_as_tikz(...).


In [13]:
print(graph.spanning_trees_as_tikz(vertex_pos, edge_styles, relabel=False))


\tikzset{edge label/.style={font=\tiny, fill=white, inner sep=1}}

\begin{tikzpicture}
	\node (e) at (0,0) {$e$};
	\node (ep) at (300:2) {$ep$};
	\node (es) at (240:2) {$es$};
	\begin{scope}[every node/.style=edge label]
		\draw[->] (ep) -- node {$k_{4}$} (e);
		\draw[->] (es) -- node {$k_{2}$} (e);
	\end{scope}
\end{tikzpicture}
%
\begin{tikzpicture}
	\node (e) at (0,0) {$e$};
	\node (ep) at (300:2) {$ep$};
	\node (es) at (240:2) {$es$};
	\begin{scope}[every node/.style=edge label]
		\draw[->] (ep) -- node {$k_{4}$} (e);
		\draw[->] (es) -- node {$k_{cat}$} (ep);
	\end{scope}
\end{tikzpicture}
%
\begin{tikzpicture}
	\node (e) at (0,0) {$e$};
	\node (ep) at (300:2) {$ep$};
	\node (es) at (240:2) {$es$};
	\begin{scope}[every node/.style=edge label]
		\draw[->] (e) to[out=330,in=90] node {$k_{3} p$} (ep);
		\draw[->] (es) -- node {$k_{2}$} (e);
	\end{scope}
\end{tikzpicture}
%
\begin{tikzpicture}
	\node (e) at (0,0) {$e$};
	\node (ep) at (300:2) {$ep$};
	\node (es) at (240:2) {$es$};
	\begin{scope}[every node/.style=edge label]
		\draw[->] (e) to[out=330,in=90] node {$k_{3} p$} (ep);
		\draw[->] (es) -- node {$k_{cat}$} (ep);
	\end{scope}
\end{tikzpicture}
%
\begin{tikzpicture}
	\node (e) at (0,0) {$e$};
	\node (ep) at (300:2) {$ep$};
	\node (es) at (240:2) {$es$};
	\begin{scope}[every node/.style=edge label]
		\draw[->] (e) to[out=210,in=90] node {$k_{1} s$} (es);
		\draw[->] (es) -- node {$k_{cat}$} (ep);
	\end{scope}
\end{tikzpicture}
%
\begin{tikzpicture}
	\node (e) at (0,0) {$e$};
	\node (ep) at (300:2) {$ep$};
	\node (es) at (240:2) {$es$};
	\begin{scope}[every node/.style=edge label]
		\draw[->] (e) to[out=210,in=90] node {$k_{1} s$} (es);
		\draw[->] (ep) -- node {$k_{4}$} (e);
	\end{scope}
\end{tikzpicture}

Basis element

graph.basis_element()

The method enzymegraph.basis_element() outputs a dictionary of the basis element of the graph.


In [14]:
for var, basis_elem in graph.basis_element().items():
    print("%s = %s" % (var, basis_elem))


ep = k_1*k_cat*s + k_2*k_3*p + k_3*k_cat*p
e = k_2*k_4 + k_4*k_cat
es = k_1*k_4*s

graph.qssa_replacements(...)

The method enzymegraph.qssa_replacements(...) outputs each the quasi-steady-state solution for each intermediate, given in terms of the total enzyme (the parameter for which should be given as the first term).


In [15]:
et = symbols('e_t', positive = True)
qssa_reps = graph.qssa_replacements(et)

for var, rep in qssa_reps.items():
    print("%s = %s" % (var, rep))


ep = e_t*(k_1*k_cat*s + k_2*k_3*p + k_3*k_cat*p)/(k_1*k_4*s + k_1*k_cat*s + k_2*k_3*p + k_2*k_4 + k_3*k_cat*p + k_4*k_cat)
e = e_t*(k_2*k_4 + k_4*k_cat)/(k_1*k_4*s + k_1*k_cat*s + k_2*k_3*p + k_2*k_4 + k_3*k_cat*p + k_4*k_cat)
es = e_t*k_1*k_4*s/(k_1*k_4*s + k_1*k_cat*s + k_2*k_3*p + k_2*k_4 + k_3*k_cat*p + k_4*k_cat)

These replacements can be used to substitute into the original ODE system.

e.g. $ \frac{dp}{dt} = k_4 ep - k_3 p e $


In [16]:
dpdt = k4*ep - k3*p*e

print("dp/dt = %s" % dpdt.subs(qssa_reps))


dp/dt = -e_t*k_3*p*(k_2*k_4 + k_4*k_cat)/(k_1*k_4*s + k_1*k_cat*s + k_2*k_3*p + k_2*k_4 + k_3*k_cat*p + k_4*k_cat) + e_t*k_4*(k_1*k_cat*s + k_2*k_3*p + k_3*k_cat*p)/(k_1*k_4*s + k_1*k_cat*s + k_2*k_3*p + k_2*k_4 + k_3*k_cat*p + k_4*k_cat)

In [16]: