means.simulation package

Submodules

Simulation Descriptors

Descriptors that are local to the simulation package

class means.simulation.descriptors.PerturbedTerm(ode_term, parameter, delta=0.01)[source]

Bases: means.core.descriptors.Descriptor

A Descriptor term that describes a particular object represents the sensitivity of some ODE term with respect to some parameter. In other words, sensitivity term describes \(s_{ij}(t) = \frac{\partial y_i(t)}{\partial p_j}\) where \(y_i\) is the ODE term described above and \(p_j\) is the parameter.

This class is used to describe sensitivity trajectories returned by means.simulation.simulate.Simulation

Parameters:
  • ode_term (ODETermBase) – the ode term whose sensitivity is being computed
  • parameter (sympy.Symbol) – parameter w.r.t. which the sensitivity is computed
delta
ode_term
parameter
class means.simulation.descriptors.SensitivityTerm(ode_term, parameter)[source]

Bases: means.core.descriptors.Descriptor

A Descriptor term that describes a particular object represents the sensitivity of some ODE term with respect to some parameter. In other words, sensitivity term describes \(s_{ij}(t) = \frac{\partial y_i(t)}{\partial p_j}\) where \(y_i\) is the ODE term described above and \(p_j\) is the parameter.

This class is used to describe sensitivity trajectories returned by means.simulation.simulate.Simulation

Parameters:
  • ode_term (ODETermBase) – the ode term whose sensitivity is being computed
  • parameter (sympy.Symbol) – parameter w.r.t. which the sensitivity is computed
mathtext()[source]
ode_term
parameter
classmethod to_yaml(dumper, data)[source]
yaml_tag = '!sensitivity-term'

Simulate

This part of the package provides utilities for simulate the dynamic of an ODEProblem.

A wide range of numerical solver are available:

>>> from means import Simulation
>>> print Simulation.supported_solvers()

In order to simulate a system, it is necessary to provide values for the initial conditions and parameters (constants):

>>> from means import mea_approximation
>>> from means.examples.sample_models import MODEL_P53
>>> from means import Simulation
>>> import numpy as np
>>>
>>> ode_problem = mea_approximation(MODEL_P53,max_order=2)
>>> # We provide initial conditions, constants and time range
>>> RATES = [90, 0.002, 1.7, 1.1, 0.93, 0.96, 0.01]
>>> INITIAL_CONDITIONS = [70, 30, 60]
>>> TMAX = 40
>>> TIME_RANGE = np.arange(0, TMAX, .1)
>>> #This is where we simulate the system to obtain trajectories
>>> simulator = Simulation(ode_problem)
>>> trajectories = simulator.simulate_system(RATES, INITIAL_CONDITIONS, TIME_RANGE)

A TrajectoryCollection object (see trajectory) is created. See the documentation of Simulation for additional information.


class means.simulation.simulate.Simulation(problem, solver='ode15s', **solver_options)[source]

Bases: means.io.serialise.SerialisableObject

Class that allows to perform simulations of the trajectories for a particular problem. Implements all ODE solvers supported by Assimulo package.

Parameters:
problem
simulate_system(parameters, initial_conditions, timepoints)[source]

Simulates the system for each of the timepoints, starting at initial_constants and initial_values values

Parameters:
  • parameters – list of the initial values for the constants in the model. Must be in the same order as in the model
  • initial_conditions – List of the initial values for the equations in the problem. Must be in the same order as these equations occur. If not all values specified, the remaining ones will be assumed to be 0.
  • timepoints – A list of time points to simulate the system for
Returns:

a list of Trajectory objects, one for each of the equations in the problem

Return type:

list[Trajectory]

solver
solver_options
classmethod supported_solvers()[source]

List the supported solvers for the simulations.

>>> Simulation.supported_solvers()
['cvode', 'dopri5', 'euler', 'lsodar', 'radau5', 'rodas', 'rungekutta34', 'rungekutta4', 'ode15s']
Returns:the names of the solvers supported for simulations
classmethod to_yaml(dumper, data)[source]
yaml_tag = '!simulation'
class means.simulation.simulate.SimulationWithSensitivities(problem, solver='ode15s', **solver_options)[source]

Bases: means.simulation.simulate.Simulation

A similar class to it’s baseclass Simulation. Performs simulations of the trajectories for each of the ODEs in given problem and performs sensitivity simulations for the problem’s parameters.

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • compute_sensitivities – Whether the model should test parameter sensitivity or not
  • solver (basestring) –

    the solver to use. Currently, the solvers that available are:

    ‘cvode’
    sundials CVode solver, as implemented in assimulo.solvers.sundials.CVode
    ‘ode15s:
    sundials CVODE solver, with default parameters set to mimick the MATLAB’s See ODE15sMixin for the list of these parameters.

    The list of these solvers is always accessible at runtime from SimulationWithSensitivities.supported_solvers() method.

  • solver_options – options to set in the solver. Consult Assimulo documentation for available options for information on specific options available.
simulate_system(parameters, initial_conditions, timepoints)[source]

Simulates the system for each of the timepoints, starting at initial_constants and initial_values values

Parameters:
  • parameters – list of the initial values for the constants in the model. Must be in the same order as in the model
  • initial_conditions – List of the initial values for the equations in the problem. Must be in the same order as these equations occur. If not all values specified, the remaining ones will be assumed to be 0.
  • timepoints – A list of time points to simulate the system for
Returns:

a list of TrajectoryWithSensitivityData objects, one for each of the equations in the problem

Return type:

list[TrajectoryWithSensitivityData]

classmethod supported_solvers()[source]

List the supported solvers for the simulations.

>>> SimulationWithSensitivities.supported_solvers()
['cvode', 'ode15s']
Returns:the names of the solvers supported for simulations

Solvers

This part of the package provides wrappers around Assimulo solvers.

class means.simulation.solvers.CVodeMixin[source]

Bases: means.simulation.solvers.UniqueNameInitialisationMixin, object

classmethod unique_name()[source]
class means.simulation.solvers.CVodeSolver(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.simulation.solvers.SolverBase, means.simulation.solvers.CVodeMixin

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
class means.simulation.solvers.CVodeSolverWithSensitivities(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.simulation.solvers.SensitivitySolverBase, means.simulation.solvers.CVodeMixin

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
class means.simulation.solvers.Dopri5Solver(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.simulation.solvers.SolverBase, means.simulation.solvers.UniqueNameInitialisationMixin

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
classmethod unique_name()[source]
class means.simulation.solvers.ExplicitEulerSolver(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.simulation.solvers.SolverBase, means.simulation.solvers.UniqueNameInitialisationMixin

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
simulate(timepoints)[source]
classmethod unique_name()[source]
class means.simulation.solvers.LSODARSolver(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.simulation.solvers.SolverBase, means.simulation.solvers.UniqueNameInitialisationMixin

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
classmethod unique_name()[source]
class means.simulation.solvers.ODE15sLikeSolver(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.simulation.solvers.SolverBase, means.simulation.solvers.ODE15sMixin

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
class means.simulation.solvers.ODE15sMixin[source]

Bases: means.simulation.solvers.CVodeMixin

A CVODE solver that mimicks the parameters used in ode15s solver in MATLAB.

The different parameters that are set differently by default are:

discr
Set to 'BDF' by default
atol
Set to 1e-6
rtol
Set to 1e-3
ATOL = 1e-06
MINH = 5.684342e-14
RTOL = 0.001
classmethod unique_name()[source]
class means.simulation.solvers.ODE15sSolverWithSensitivities(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.simulation.solvers.SensitivitySolverBase, means.simulation.solvers.ODE15sMixin

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
class means.simulation.solvers.Radau5Solver(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.simulation.solvers.SolverBase, means.simulation.solvers.UniqueNameInitialisationMixin

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
classmethod unique_name()[source]
class means.simulation.solvers.RodasSolver(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.simulation.solvers.SolverBase, means.simulation.solvers.UniqueNameInitialisationMixin

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
classmethod unique_name()[source]
class means.simulation.solvers.RungeKutta34Solver(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.simulation.solvers.SolverBase, means.simulation.solvers.UniqueNameInitialisationMixin

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
classmethod unique_name()[source]
class means.simulation.solvers.RungeKutta4Solver(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.simulation.solvers.SolverBase, means.simulation.solvers.UniqueNameInitialisationMixin

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
simulate(timepoints)[source]
classmethod unique_name()[source]
class means.simulation.solvers.SensitivitySolverBase(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.simulation.solvers.SolverBase

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
class means.simulation.solvers.SolverBase(problem, parameters, initial_conditions, starting_time=0.0, **options)[source]

Bases: means.util.memoisation.MemoisableObject

This acts as a base class for ODE solvers used in means. It wraps around the solvers available in :module:`assimulo` package, and provides some basic functionality that allows solvers be used with means objects.

Parameters:
  • problem (ODEProblem) – Problem to simulate
  • parameters (iterable) – Parameters of the solver. One entry for each constant in problem
  • initial_conditions (iterable) – Initial conditions of the system. One for each of the equations. Assumed to be zero, if not specified
  • starting_time (float) – Starting time for the solver, defaults to 0.0
  • options – Options to be passed to the specific instance of the solver.
simulate(timepoints)[source]

Simulate initialised solver for the specified timepoints

Parameters:timepoints – timepoints that will be returned from simulation
Returns:a list of trajectories for each of the equations in the problem.
exception means.simulation.solvers.SolverException(message, base_exception=None)[source]

Bases: exceptions.Exception

base_exception
class means.simulation.solvers.UniqueNameInitialisationMixin[source]

Bases: object

classmethod unique_name()[source]
means.simulation.solvers.available_solvers(with_sensitivity_support=False)[source]
means.simulation.solvers.parse_flag(exception_message)[source]

Parse the flag from the solver exception. e.g.

>>> parse_flag("Exception: Dopri5 failed with flag -3")
-3
Parameters:exception_message (str) – message from the exception
Returns:flag id
Return type:int

Gillespie Stochastic Simulation Algorithm

This part of the package provides a simple implementation of GSSA. This is designed for experimental purposes much more than for performance. If you would like to use SSA for parameter inference, or for high number of species, there are many superior implementations available.

class means.simulation.ssa.SSASimulation(stochastic_problem, n_simulations, random_seed=None)[source]

Bases: means.io.serialise.SerialisableObject

A class providing an implementation of the exact Gillespie Stochastic Simulation Algorithm [Gillespie77].

>>> from means.examples import MODEL_P53
>>> from means import StochasticProblem, SSASimulation
>>> import numpy as np
>>> PROBLEM = StochasticProblem(MODEL_P53)
>>> RATES = [90, 0.002, 1.7, 1.1, 0.93, 0.96, 0.01]
>>> INITIAL_CONDITIONS = [70, 30, 60]
>>> TIME_RANGE = np.arange(0, 40, .1)
>>> N_SSA = 10
>>> ssas = SSASimulation(PROBLEM, N_SSA)
>>> mean_trajectories = ssas.simulate_system(RATES, INITIAL_CONDITIONS, TIME_RANGE)
Parameters:
  • stochastic_problem
  • n_simulations
  • random_seed
simulate_system(parameters, initial_conditions, timepoints, max_moment_order=1, number_of_processes=1)[source]

Perform Gillespie SSA simulations and returns trajectories for of each species. Each trajectory is interpolated at the given time points. By default, the average amounts of species for all simulations is returned.

Parameters:
  • parameters – list of the initial values for the constants in the model. Must be in the same order as in the model
  • initial_conditions – List of the initial values for the equations in the problem. Must be in the same order as these equations occur.
  • timepoints – A list of time points to simulate the system for
  • number_of_processes – the number of parallel process to be run
  • max_moment_order – the highest moment order to calculate the trajectories to. if set to zero, the individual trajectories will be returned, instead of the averaged moments.

E.g. a value of one will return means, a values of two, means, variances and covariance and so on.

Returns:a list of Trajectory one per species in the problem, or a list of lists of trajectories (one per simulation) if return_average == False.
Return type:list[Trajectory]
means.simulation.ssa.multiprocessing_apply_ssa(x)[source]

Used in the SSASimulation class. Needs to be in global scope for multiprocessing module to pick it up

means.simulation.ssa.multiprocessing_pool_initialiser(population_rates_as_function, change, species, initial_conditions, t_max, seed)[source]

Trajectories

This part of the package provide convenient utilities to manage trajectories. A Trajectory object is generally a time series containing the values of a given moment (e.g. mean, variance, ...) over a time range. Trajectories are typically returned by simulations (see simulate and ssa), or from observation/measurement.

The TrajectoryCollection class is a container of trajectories. It can be used like other containers such as lists.

Both ~means.simulation.trajectory.TrajectoryCollection and ~means.simulation.trajectory.Trajectory have there own .plot() method to help representation.

class means.simulation.trajectory.Trajectory(timepoints, values, description)[source]

Bases: means.io.serialise.SerialisableObject

A single simulated or observed trajectory for an ODE term.

Parameters:
  • timepoints (iterable) – timepoints the trajectory was simulated for
  • values (iterable) – values of the curve at each of the timepoints
  • description (Descriptor) – description of the trajectory
description

Description of this trajectory. The same description as the description for particular ODE term.

Return type:Descriptor
plot(*args, **kwargs)[source]

Plots the trajectory using matplotlib.pyplot.

Parameters:
  • args – arguments to pass to plot()
  • kwargs – keyword arguments to pass to plot()
Returns:

the result of the matplotlib.pyplot.plot() function.

png
resample(new_timepoints, extrapolate=False)[source]

Use linear interpolation to resample trajectory values. The new values are interpolated for the provided time points. This is generally before comparing or averaging trajectories.

Parameters:
  • new_timepoints – the new time points
  • extrapolate – whether extrapolation should be performed when some new time points are out of the current time range. if extrapolate=False, it would raise an exception.
Returns:

a new trajectory.

Return type:

Trajectory

set_description(description)[source]
svg
timepoints

The timepoints trajectory was simulated for.

Return type:numpy.ndarray
to_csv(file)[source]

Write this trajectory to a csv file with the headers ‘time’ and ‘value’.

Parameters:file (file) – a file object to write to
Returns:
classmethod to_yaml(dumper, data)[source]
values

The values for each of the timepoints in timepoints.

Return type:numpy.ndarray
yaml_tag = u'!trajectory'
class means.simulation.trajectory.TrajectoryCollection(trajectories)[source]

Bases: means.io.serialise.SerialisableObject

A container of trajectories with representation functions for matplotlib and IPythonNoteBook. In most cases, it simply behaves as list.

plot(legend=True)[source]
png
svg
to_csv(file)[source]

Write all the trajectories of a collection to a csv file with the headers ‘description’, ‘time’ and ‘value’.

Parameters:file (file) – a file object to write to
Returns:
classmethod to_yaml(dumper, data)[source]
trajectories

Return a list of all trajectories in the collection :rtype: list[Trajectory]

yaml_tag = '!trajectory-collection'
class means.simulation.trajectory.TrajectoryWithSensitivityData(timepoints, values, description, sensitivity_data)[source]

Bases: means.simulation.trajectory.Trajectory

An extension to Trajectory that provides data about the sensitivity of said trajectory as well.

Parameters:
  • timepoints (numpy.ndarray) – timepoints the trajectory was simulated for
  • values (numpy.ndarray) – values of the curve at each of the timepoints
  • description (Descriptor) – description of the trajectory
  • sensitivity_data (list[Trajectory]) – a list of Trajectory objects signifying the sensitivity change over time for each of the parameters.
classmethod from_trajectory(trajectory, sensitivity_data)[source]
plot_perturbations(parameter, delta=0.0001, *args, **kwargs)[source]
sensitivity_data

THe sensitivity data for the trajectory

Return type:list[Trajectory]
classmethod to_yaml(dumper, data)[source]
yaml_tag = '!trajectory-with-sensitivity'
means.simulation.trajectory.perturbed_trajectory(trajectory, sensitivity_trajectory, delta=0.0001)[source]

Slightly perturb trajectory wrt the parameter specified in sensitivity_trajectory.

Parameters:
  • trajectory (Trajectory) – the actual trajectory for an ODE term
  • sensitivity_trajectory (Trajectory) – sensitivity trajectory (dy/dpi for all timepoints t)
  • delta (float) – the perturbation size
Returns:

Trajectory

Module contents

Routines for stochastic and deterministic simulation.