The pybrn module reference

brn.brn

provides the class Brn, representing a biochemical reaction network.

class brn.brn.Brn(stoich, species=None, reactions=None, parameters=[], variables=[], rates={}, ics={}, parvals={}, varvals={}, scale=1.0)

This class represents a biochemical reaction network, specifically the ODE of the form xdot = N * v(x), where N is the stoichiometric matrix, v the reaction rate vector, and x the species concentration vector.

methods: Brn.addelements - add reactions, species, and/or parameters to the BRN Brn.conservation_analysis - analyse conservation relations Brn.eval - evaluate a symbolic expressions in the context of this BRN Brn.jacobian - compute a Jacobian matrix Brn.rate_fun - generate a lambda function for the rate vector Brn.scale_fun - generate a lambda function for the scale matrix Brn.set - set values for this BRN Brn.showreact - display a reaction as string Brn.simulate - simulate the network’s ODE Brn.steadystate - compute a steady state Brn.steadystate_branch - compute a branch of steady states by varying

one parameter

properties: Brn.has_conservation - flag indicating whether conserved moieties are present Brn.has_scale - flag indicating whether a non-identiy scale factor is present

addelements(reactions=[], stoich={}, species=[], parameters=[], values={}, scale=1.0)

add the given reactions, species, and parameters (sequences of strings) to the network reaction rates and other values may be given in the dict ‘’values’‘

If reactions is non-empty, then ‘’stoich’’ is a dictionary of dictionaries describing the stoichiometry of the reactions. For example, if reaction ‘’v’’ is added, and the network contains species A which is converted to B through v, then stoich = {‘v’:{‘A’:-1, ‘B’:1}}

conservation_analysis(store=True)

do a conservation analysis of the network.

generates the attribute ‘conserv’ in this Brn object which is a brn.conserv.Conservation object, containing all information about conservation relations. The generated Conservation object is also returned by the method.

See the documentation of brn.conserv.Conservation for further information.

eval(expr, valuedict={}, p=None, x=None, values={}, level='all')

evaluate a symbolic expression, given as string ‘expr’, for this BRN state and parameter values may be given in the dictionary ‘values’, or

in the vectors ‘x’ for state variables and ‘p’ for parameters. the values in ‘values’ take precedence over the vectors. for unspecified variables default values are used.

see the documentation of brneval.ExpressionEvaluator for more info.

evalparam(expr, p=None, values={})
evaluate a symbolic expression, given as string ‘expr’, which only
depends on the network’s parameters, for this BRN
parameter values may be given in the dictionary ‘values’ or in the vector ‘p’.
for unspecified variables default values are used.
has_conservation

returns true if the network has conserved moieties, false otherwise.

has_scale

returns a Boolean flag which is true if the network has defined a non-identity scale factor (xdot = scale * N * v(x)), and false otherwise.

jacobian(fun, var, returntype='fun')

compute a jacobian of the brn, either as lambda function or as symbolic expression

arguments: fun - symbolic expression which the Jacobian should be taken of var - variables in which the differentiation is done returntype - specify one of ‘fun’ or ‘str’

Usage examples: >>> jac = net.jacobian(net.reactions,net.species) >>> jac(net.eval(net.species),net.eval(net.parameters))

linearapproximation(x=None, p=None, values={}, inputs=None, outputs=None)

compute matrices for a linear approximation this Brn.

Arguments: values - dict to specify custom species and parameter values x, p - custom species, parameter values, takes precedence over values inputs - list of parameters to be considered as input variables outputs - list of symbolic expressions to be considered as output variables

Compute the linear approximations A = df/dz, B = df/du, and C = dh/dx, where f is the network’s ODE right hand side, z the state variable, u the input variables, and h the output functions, evaluated at the specified species concentrations and parameter values

Returns (A,B,C), where B and / or C are None if inputs / outputs is.

rate_fun(values={}, returntype='fun')

generate a lambda function for the reaction rates v(x,p) argument to the generated function are lists of species and parameter values

Custom mathematical expressions for reaction rate or variables can be passed in ‘values’.

Options for ‘returntype’: fun - return a callable function str - return the string of the function definition

usage: >>> fun = net.rate_fun() >>> fun([1],[]) -1.0

scale_fun()

generate a lambda function for the scaling factor argument to the generated function is a list of parameter values

usage: >>> fun = net.scale_fun() >>> fun([]) 1

set(values, strict=True)

Assign values to symbols in the network from the dict values

Usage example: >>> net.set({‘k1’:1}) # sets parameter k1 to 1

Set strict to False to avoid exceptions when key in values are not known in the network.

showreact(react, printstr=True)

generate nice string representations of the reaction ‘react’ if ‘printstr’ is True, output the generated strings via the print method

usage: >>> strings = net.showreact(‘v1’)

simulate(t, p=None, x0=None, values={}, method=<class 'brn.simulation.VodeSimulator'>, **args)

simulate the network

arguments: t - list of time points for which simulation data should be generated method - which implementation of brn.simulation.Simulator to use args - keyword arguments which are passed to the simulator upon initialization

see e.g. help(brn.simulation.Simulator) and subclasses
steadystate(p=None, x0=None, values={})

try to compute a steady state of the Brn

output: xs - result of steady state computation

custom initial conditions for species and parameter values can be given in ‘x0’ or ‘p’, respectively state and parameter values may also be updated via the dict ‘values’ information from ‘x0’ and ‘p’ takes precedence over ‘values’

steadystate_branch(par, start, stop, step=0, p=None, x0=None, values={})

simplistic computation of a steady state branch: varies the value of parameter ‘par’ from ‘start’ to ‘stop’, increasing the value in steps of length ‘step’ if ‘step’ is 0, make a total of 100 steps from start to stop

custom initial conditions for species and parameter values can be given in ‘x0’ or ‘p’, respectively state and parameter values may also be updated via the dict ‘values’ information from ‘x0’ and ‘p’ takes precedence over ‘values’

output: plist - list of parameter values for which steady states have been computed xslist - list of steady states, corresponding to parameter values in plist

brn.brn.make_massaction_network(stoich, species=None, reactions=None, parameters=None, paramvalues=None, ics=None)

generate a Brn with standard mass action reaction rates from the stoichiometric matrix ‘stoich’.

This function builts the network with irreversible reactions only!

Optional arguments can be given as list of species names, reaction names, parameter names, numerical parameter values, and numerical initial conditions for species concentrations. If species, reaction, or parameter names are not given, they default to x1, x2, ..., v1, v2, ..., and k1, k2, ..., respectively. If parameter values or initial conditions are not given, they default to None.

Note that parameter values and initial conditions can be changed later via the Brn.set() method.

brn.brncomp

brncomp - module for symbolic and derivative computations in biochemical networks

Methods: brncomp.jacobian - computation of Jacobians using sympy brncomp.linearapproximation - compute matrices for a linear approximation

brn.brncomp.jacobian(f, x, vars=[], valuedict={}, returntype='fun')

compute the Jacobian of f with respect to x f and x are string lists vars is a list of lists of independent variables, which should be maintained as variables in the resulting Jacobian. if vars is empty, take vars=[x] valuedict is a dict of symbols to be substituted by either a number or

a symbolic expression, with a substitution depth up to two, i.e. values in valuedict may themselves contain keys in valuedict, but these keys most contain only numbers or symbols from x or param.
returntype may be one of:
“fun” - return a lambda function for the jacobian, arguments are lists of independent
variables as given in the tuple ‘vars’

“str” - return the jacobian as list of symbolic expressions “sym” - return as sympy expressions

usage example: >>> jacobian([‘f’],[‘x’],[[‘x’],[‘k’]],valuedict={‘f’:’k*x**2’},returntype=”str”) [[‘2*k*x’]] >>> fun = jacobian([‘f’],[‘x’],[[‘x’],[‘k’]],valuedict={‘f’:’k*b’,’b’:’x**2’},returntype=”fun”) >>> fun([2],[2]) [[8]]

brn.brncomp.linearapproximation(net, x=None, p=None, inputs=None, outputs=None, returntype='fun')

compute matrices for a linear approximation of Brn ‘net’.

Arguments: returntype - one of

‘fun’: resulting matrix M is a function to evaluate as M(x,p) ‘num’: resulting matrix M is a numpy array, evaluated at species ‘x’ and parameters ‘p’

using default values for ‘x’ or ‘p’ if they are None

x, p - custom species, parameter values inputs - list of parameters to be considered as input variables outputs - list of symbolic expressions to be considered as output variables

Compute the linear approximations A = df/dz, B = df/du, and C = dh/dx, where f is the network’s ODE right hand side, z the state variable, u the input variables, and h the output functions.

Returns (A,B,C), where B and / or C are None if inputs / outputs is.

brn.brncomp.variable_substitution(vardict, params)

Substitute variable references in symbolic expressions. ‘vardict’ is a dict with variable names as keys and formulas as values. Formulas may depend on other variables and symbolic parameters given in the list ‘params’. Tries to substitute all references to other variables in the variable formulas by the corresponding expressions.

May raise ValueError if problems arise, e.g. due to an algebraic loop.

brn.brneval

numerical evaluation of Brns

class brn.brneval.BrnEvaluator(net, p=None, x0=None, values={})

Abstract class for numerical computation tasks: - interface for evaluating the network’s ODE right hand side and Jacobian - useful for simulation and steady state computation

Classes that do numerical computation can inherit from this class.

convert_points(z)

convert the list of vectors ‘z’ in the internal, reduced coordinates, to the original coordinate system

returns: x - a list of vectors corresponding to z in the original coordinates

ode_jac(z)

evaluate the Jacobian of the BRN’s ODE right hand side This method must be implemented by classes derived from the BrnEvaluator class.

ode_rhs(z)

evaluate the BRN’s ODE right hand side This method must be implemented by classes derived from the BrnEvaluator class.

set(p=None, x0=None, values={})

sets parameter values to ‘p’ and initial point to ‘x0’ for this brn evaluation problem numerical values for x and p can also be given in dict values, these override settings from ‘p’ and ‘x0’

class brn.brneval.CythonBrnEvaluator(net, p=None, x0=None, values={}, precompiled=None)

Subclass of BrnEvaluator that generates a Cython-compiled module to evaluate the network’s right hand side. This should be faster than the PythonBrnEvaluator for many subsequent evaluations of the same network.

Requires that Cython ( http://cython.org ) is installed.

compile(precompiled=None)

Compile the network into a loadable module with cython.

If the same network has been compiled previously, give the full path (without extension) to the resulting module in the argument ‘precompiled’, the compilation is then skipped.

Returns the full path to the resulting module (without extension).

ode_jac(z)

call Cython implementation of the network’s rhs Jacobian

ode_rhs(z)

call Cython implementation of the network’s rhs fun = S N v(x,p)

set(p=None, x0=None, values={})

sets parameter values to ‘p’ and initial point to ‘x0’ for this brn evaluation problem

class brn.brneval.ExpressionEvaluator(net, p=None, x=None, values={}, typestrict=True, level='all')

Class to evaluate arbitrary expressions in the context of a given BRN.

evalexpression(expr)

evaluate a symbolic expression, given as string ‘expr’, with this evaluator. may also be a sequence of symbolic expressions

set(p=None, x=None, values={})

update internal values of this evaluator to values given in vectors p, x0, and dict values. values takes precedence over p, x0.

Returns the updated ExpressionEvaluator object.

class brn.brneval.PythonBrnEvaluator(net, p=None, x0=None, values={})

Subclass of BrnEvaluator implementing the rhs functions directly in Python. A simple but slow implementation without any further dependencies.

ode_jac(z)

Python implementation of the network’s rhs Jacobian jacfun = d ode_rhs / dx

ode_rhs(z)

Python implementation of the network’s rhs fun = S N v(x,p)

brn.conserv

conserv - provides methods for the analysis of conservation relations

class brn.conserv.Conservation(net)

provide tools for the analysis of networks with conserved moieties

methods: full2reduced_fun - generate a lambda function for conversion from full to reduced coordinates reduced2full_fun - generate a lambda function for conversion from reduced to full coordinates

full2reduced_fun()

return a lambda function for the conversion from the original coordinate system to the reduced coordinates

the formula is z = L’ x

reduced2full_fun(x0)

return a lambda function for the conversion from the reduced coordinate system to the original coordinates

the formula is x = L z + G G’ x0

brn.mathutils

mathutils module: provides basic mathematics functionality for pybrn

functions: rank - compute the matrix rank

brn.mathutils.rank(A)

compute the rank of matrix A

brn.sbmlimport

methods to import sbml to pybrn data structure

class brn.sbmlimport.FunctionHelper(sbml_model)

Helper class to deal with SBML function definitions

process(expr)

Process expr by replacing all function calls with the corresponding formula

brn.sbmlimport.fromSBML(filename, duplicatelocalparameters=False)

Import an SBML file as BRN. Usage: net = fromSBML(filename) where filename is the SBML file and net is the resulting BRN.

Set duplicatelocalparameters=True if local parameters with equal names should NOT be treated as the same parameters.

brn.sbmlimport.fromSBMLString(sbmlstring, duplicatelocalparameters=False)

Generate a BRN from an SBML string Usage: net = fromSBML(sbmlstring) where sbmlstring is a string in SBML format and net is the resulting BRN.

Set duplicatelocalparameters=True if local parameters with equal names should NOT be treated as the same parameters.

brn.sbmlimport.readSBML(document, duplicatelocalparameters=False, undefined=nan)

Convert an SBML document to a BRN. Usage: net = fromSBML(document) where document is a libsbml.SBMLDocument object and net is the resulting BRN.

Set duplicatelocalparameters=True if local parameters with equal names should NOT be treated as the same parameters.

‘undefined’ is used for any quantities without an explicit assignment in the SBML file, apart from compartments, which use 1.0 in this case.

brn.steadystate

modul for computation of steady states in biochemical networks

class brn.steadystate.SteadyState(net, p=None, x0=None, values={}, evaluator=<class 'brn.brneval.PythonBrnEvaluator'>)

steady state computation for biochemical networks

important methods: - SteadyState.set: set values for parameters and/or initial conditions - SteadyState.compute: solve for the steady state

compute()

try to solve the steady state computation problem using scipy.optimize.fsolve returns xs, info, with: xs - result of computation info - dictionary of additional information, see scipy.optimize.fsolve documentation

brn.steadystate.steadystate_branch(net, par, start, stop, step, p, x0)

compute a steady state branch for network ‘net’, with parameter ‘par’ varying from ‘start’ to ‘stop’. the parameter is varied in steps of length ‘step’ additional arguments: p - parameter values to be used x0 - initial point for the computation

output: plist - list of values for parameter ‘par’ xslist - list of steady states corresponding to the parameter values in ‘plist’ info - last output of the call to SteadyState.compute()

brn.simulation

modul for numerical simulation of biochemical networks

class brn.simulation.GillespieSimulator(net, times=[0, 1], x0=None, p=None, values={}, evaluator=<class 'brn.brneval.PythonBrnEvaluator'>, outputs=('species',))

implements network simulation with the Gillespie algorithm

compute_outputs(xlist, outputs=('species', ))

from the list of species concentrations ‘xlist’, compute variables according to the sequence ‘outputs’. Elements of ‘outputs’ may be ‘species’, ‘reactions’, ‘variables’, or named network elements. The outputs are returned as tuple in the order as they are given in ‘outputs’.

run(times=None, outputs=None)

Method to solve the simulation task. Simulation data is generated at the time points given in the list ‘times’, and at default points if times=None.

outputs - BRN variables to compute during simulation

output: t - list of time points from simulation further outputs as specified in ‘outputs’, default is list of state values representing the trajectory

class brn.simulation.Simulator(net, times=[0, 1], x0=None, p=None, values={}, evaluator=<class 'brn.brneval.PythonBrnEvaluator'>, outputs=('species',))

abstract class for simulation of biochemical networks

Simulation classes should inherit from this class and at least redefine the run() method.

compute_outputs(z, outputs=('species', ))

from the list of vectors ‘z’ in internal coordinates, compute variables according to the sequence ‘outputs’. Elements of ‘outputs’ may be ‘species’, ‘reactions’, ‘variables’, or named network elements. The outputs are returned as tuple in the order as they are given in ‘outputs’.

convert_trajectory(z)

convert the list of vectors ‘z’ in the internal, reduced coordinates, to the original coordinate system

returns: x - a list of vectors corresponding to z in the original coordinates

initialize(x0=None, p=None, values={})

set the initial condition for this simulation task to ‘x0’, a list or array of initial condition values. optionally also set the parameter values to ‘p’

run()

Method to solve the simulation task. This method must be implemented by classes derived from the Simulator class.

class brn.simulation.SundialsSimulator(net, times=[0, 1], x0=None, p=None, values={}, precompiled=None, filename=None, outputs=('species',), library_dirs=[], noload=False)

implements network simulation via cython and the sundials C library

compile(precompiled=None, filename=None, library_dirs=[])

Compile a sundials wrapper for this network.

If the same network has been compiled previously, give the full path (without extension) to the resulting module in the argument ‘precompiled’, the compilation is then skipped.

Returns the full path to the resulting module (without extension).

do_compilation(fname, library_dirs=[])

compile a sundials wrapper for this network and store in ‘filename’ (platform dependent suffix is added automatically)

intermediate build files are stored in the same directory as filename

initialize(x0=None, p=None, values={})

set the initial condition for this simulation task to ‘x0’, a list or array of initial condition values. optionally also set the parameter values to ‘p’

run(times=None, outputs=None)

Method to solve the simulation task. Simulation data is generated at the time points given in the list ‘times’, and at default points if times=None.

outputs - BRN variables to compute during simulation

output: t - list of time points from simulation further outputs as specified in ‘outputs’, default is list of state values representing the trajectory

class brn.simulation.VodeSimulator(net, times=[0, 1], x0=None, p=None, values={}, evaluator=<class 'brn.brneval.PythonBrnEvaluator'>, outputs=('species',))

implements network simulation via scipy.integrate

run(times=None, outputs=None)

Method to solve the simulation task. Simulation data is generated at the time points given in the list ‘times’, and at default points if times=None.

outputs - BRN variables to compute during simulation

output: t - list of time points from simulation further outputs as specified in ‘outputs’, default is list of state values representing the trajectory

class brn.simulation.VodeSimulatorStop(net, times=[0, 1], x0=None, p=None, values={}, evaluator=<class 'brn.brneval.PythonBrnEvaluator'>, outputs=('species',), stopcondition=False)

Implements network simulation via scipy.integrate. Allows the specification of stop conditions on the simulation. Simulation stops at the first time point where a simulation result has been requested and the stop condition evaluates to True.

run(times=None, outputs=None)

Method to solve the simulation task. Simulation data is generated at the time points given in the list ‘times’, and at default points if times=None.

outputs - BRN variables to compute during simulation

output: t - list of time points from simulation further outputs as specified in ‘outputs’, default is list of state values representing the trajectory

symbolic_stopcondition(x)

Returns true if simulation should be stopped for current state x, False otherwise.

brn.simulation.compile_model(net, modulename, builddir='.', library_dirs=[])

compile model ‘net’ to a shared library module for use with SundialsSimulator

modulename - name of the shared library (potentially including suffix .so) builddir - directory where build files are placed in (default is current directory)

Example: >>> compile_model(net, “net_module”) >>> sim = SundialsSimulator(net, times=np.arange(10.0), precompiled=”net_module”)

Table Of Contents

Previous topic

The pybrn Tutorial

This Page