DiffEqBiological.add_hc_template!Method
add_hc_template(reaction_network)

Adds another hc template to the network. Having more HC templates might increase run time a little bit, but might give some additional certainity that all solutions are found (this should generally not be a problem anyway)

args

  • reaction_network: a reaction network.
DiffEqBiological.add_scale_noise_param!Method
add_scale_noise_param!(network, scale_noise_name::String)

Given an AbstractReaction network, add the parameter with the passed in string as its name to the network (if it is not already defined), and register it as the noise scaling coefficient.

DiffEqBiological.add_scale_noise_param!Method
add_scale_noise_param!(network, scale_noise::Symbol)

Given an AbstractReaction network, add the parameter corresponding to the passed in symbol to the network (if it is not already defined), and register it as the noise scaling coefficient.

DiffEqBiological.addjumps!Method
addjumps!(network; build_jumps=true, build_regular_jumps=true, minimal_jumps=false)

Extend an AbstractReactionNetwork generated with the @min_reaction_network or @empty_reaction_network macros with everything needed to use jump SSA solvers.

Optional kwargs can be used to disable the construction of additional jump solver components.

Keyword arguments:

  • build_jumps: if true jump rates and affects will be calculated for use in DiffEqJump SSAs.

  • build_regular_jumps: if true a RegularJump representation of the stochastic chemical kinetics model will be calculated for use in τ-leaping methods.

  • minimal_jumps: if true ConstantRate jumps are only constructed for non-mass action jumps. (Note, mass action jumps are still resolved within any jump simulation. This option simply speeds up the construction of the jump problem since it avoids building redundant ConstantRate jumps that encode MassActionJumps, which are subsequently ignored within jump simulations.)

DiffEqBiological.addodes!Method
addodes!(network; build_jac=true, build_symfuncs=true)

Extend an AbstractReactionNetwork generated with the @min_reaction_network or @empty_reaction_network macros with everything needed to use ODE solvers.

Optional kwargs can be used to disable the construction of additional ODE solver components.

DiffEqBiological.addparam!Method
addparam!(network, paramname::String)

Given an AbstractReaction network, add the parameter with name given by the passed in string to the network (if it is not already defined).

DiffEqBiological.addparam!Method
addparam!(network, param::Symbol)

Given an AbstractReaction network, add the parameter corresponding to the passed in symbol to the network (if it is not already defined).

DiffEqBiological.addreaction!Method
addreaction!(network, rateex::Union{Expr,Symbol,Int,Float64}, rxexpr::Expr)

Given an AbstractReaction network, add a reaction with the passed in rate and reaction expressions. i.e. a reaction of the form

k*X, 2X + Y --> 2W

would have rateex=:(k*X) and rxexpr=:(2X + Y --> W),

10.5, 0 --> X

would have rateex=10.5 and rxexpr=:(0 --> X), and

k, X+X --> Z

would have rateex=:k and rxexpr=:(X+X --> Z). All normal DSL reaction definition notation should be supported.

DiffEqBiological.addreaction!Method
addreaction!(network, rateex::Union{Expr,Symbol,Int,Float64}, substrates, products)

Given an AbstractReaction network, add a reaction with the passed in rate, rateex, substrate stoichiometry, and product stoichiometry. Stoichiometries are represented as tuples of Pair{Symbol,Int}. i.e. a reaction of the form

k*X, 2X + Y --> 2W

would have rateex=:(k*X), substrates=(:X=>2, :Y=>2)andproducts=(W=>2,)`,

10.5, 0 --> X

would have rateex=10.5, substrates=() and products=(:X=>1,), and

k, X+X --> Z

would have rateex=:k, substrates=(:X=>2,) and products=(:Z=>2,). All normal DSL reaction definition notation should be supported for the rateex.

DiffEqBiological.addsdes!Method
addsdes!(network; build_jac=true, build_symfuncs=true)

Extend an AbstractReactionNetwork generated with the @min_reaction_network or @empty_reaction_network macros with everything needed to use SDE solvers.

Optional kwargs can be used to disable the construction of additional SDE solver components.

DiffEqBiological.addspecies!Method
addspecies!(network, speciesname::String)

Given an AbstractReaction network, add the species with name given by the passed in string to the network (if it is not already defined.

DiffEqBiological.addspecies!Method
addspecies!(network, speciessym::Symbol)

Given an AbstractReaction network, add the species corresponding to the passed in symbol to the network (if it is not already defined).

DiffEqBiological.bifurcation_gridMethod
bifurcation_grid(reaction_network, [solver], parameter_values, parameter_symbol, parameter_range; kwargs...)

Solves the steady states of a system for every parameter values in a given range of values. Returns a bifurcation grid type object (which can be plotted).

args

  • reaction_network: a reaction network.
  • solver (optional): a subtype of AbstractSteadyStateSolver that specifies how to solve the bifurcation diagram. Default is HCSteadyStateSolver (currently only solver avaiable).
  • parameter_values: a vector which specifies the parameter values for which the steady states are to be found. Not required in case the reaction network lacks steady states.
  • parameter_symbol: the parameter which is to be varried.
  • parameter_range: a range with all the values of the parameter for which steady states are to be solved.

kwargs

  • potential arguments for steady state solvers (however, only current solver does not have any such arguments).
DiffEqBiological.bifurcation_grid_2dMethod
bifurcation_grid_2d(reaction_network, [solver], parameter_values, parameter_symbol1, parameter_range1, parameter_symbol2, parameter_range2; kwargs...)

Ad bifurcation grid, but varries two different parameters. Returning a 2d grid with steady states in each grid point.

args

  • reaction_network: a reaction network.
  • solver (optional): a subtype of AbstractSteadyStateSolver that specifies how to solve the bifurcation diagram. Default is HCSteadyStateSolver (currently only solver avaiable).
  • parameter_values: a vector which specifies the parameter values for which the steady states are to be found. Not required in case the reaction network lacks steady states.
  • parameter_symbo1l: the first parameter which is to be varried.
  • parameter_range1: a range with all the values of the first parameter for which steady states are to be solved.
  • parameter_symbol2: the second parameter which is to be varried.
  • parameter_range2: a range with all the values of the second parameter for which steady states are to be solved.

kwargs

  • potential arguments for steady state solvers (however, only current solver does not have any such arguments).
DiffEqBiological.bifurcation_grid_diagramMethod
bifurcation_grid_diagram(reaction_network, [solver], parameter_values, parameter_symbol1, parameter_range1, parameter_symbol2, parameter_range2; kwargs...)

Varries a first parameter over a range of discrete values, for each values makes a bifurcation diagram over a second, continious range. Allows to vsiualise changes ins steady states over two parameters.

args

  • reaction_network: a reaction network.
  • solver (optional): a subtype of AbstractBifurcationSolver that specifies how to solve the bifurcation diagram. Default is HCBifurcationSolver.
  • parameter_values: a vector which specifies the point in parameter space around which the bifurcation diagram is evaluated.
  • parameter_symbol1: the first parameter which is varied, given as a symbol.
  • parameter_range1: the range over which the first parameter is varied (this is a range of discrete numbers).
  • parameter_symbol2: the second parameter which is varied, given as a symbol.
  • parameter_rang21: the range over which the second parameter is varied (this is a tuple of two numbers).

kwargs

  • stepsize = (parameterrange[2] - parameterrange[1])/200: The maximum distance between two parameter values for which steady states are calculated. Smaller values will give a smoother diagram.
  • Δp = (parameterrange[2] - parameterrange[1])/200: The distance to jump after finding a bifurcation. After discovering a bifurcation, the solver jumps ahead a distance Δp and starts back-tracking. If this distance is too small, the method may error or cause visible artifiacts in the bifurcation diagram. If it is too large, then you might jump over another bifurcation and it will be missed (this should also be farily obvious in a plot). Only applicable for the HCBifurcationSolver.
  • Δu = 1e-7: (alternatively reduced to a tenth of the minimum value of any initial solutions) used by algorithm to determine if two points are identical. Set this to less than the expected minimum distance between two different solutions.
DiffEqBiological.bifurcationsMethod
bifurcations(reaction_network, [solver], parameter_values, parameter_symbol, parameter_range; kwargs...)

Get a bifurcation diagram of the specified system.

args

  • reaction_network: a reaction network.
  • solver (optional): a subtype of AbstractBifurcationSolver that specifies how to solve the bifurcation diagram. Default is HCBifurcationSolver.
  • parameter_values: a vector which specifies the point in parameter space around which the bifurcation diagram is evaluated.
  • parameter_symbol: the parameter which is varied, given as a symbol.
  • parameter_range: the range over which the specified parameter is varied.

kwargs

  • stepsize = (parameterrange[2] - parameterrange[1])/200: The maximum distance between two parameter values for which steady states are calculated. Smaller values will give a smoother diagram.
  • Δp = (parameterrange[2] - parameterrange[1])/200: The distance to jump after finding a bifurcation. After discovering a bifurcation, the solver jumps ahead a distance Δp and starts back-tracking. If this distance is too small, the method may error or cause visible artifiacts in the bifurcation diagram. If it is too large, then you might jump over another bifurcation and it will be missed (this should also be farily obvious in a plot). Only applicable for the HCBifurcationSolver.
  • Δu = 1e-7: (alternatively reduced to a tenth of the minimum value of any initial solutions) used by algorithm to determine if two points are identical. Set this to less than the expected minimum distance between two different solutions.
DiffEqBiological.dependentsMethod
dependents(network, rxidx)

Given an AbstractReactionNetwork and a reaction index, rxidx, return a vector of symbols of species the reaction rate law depends on. i.e. for

k*W, 2X + 3Y --> 5Z + W

the returned vector would be [:W,:X,:Y].

Non-allocating, returns underlying field within the reaction_network.

DiffEqBiological.fix_parametersFunction
fix_parameters(reaction_network, [parameter_values]; t=get_polyvars(rn).t, kwargs...)

Fixes certain parameters of a reaction network, preparing it for homotopy continuation.

In the polynomial needed to be solved by homotopy continuation to find steady states solutions, both varriables and parameters are initially represented as polynomial varriables. If a parameter exists in an exponent (e.g. X^n), that parameter might be needed to be fixed to a desired value before solving. Fixing is not required, but if you solve a polynomial several times, it might increase performance.

args

  • reaction_network: a reaction network.
  • parameter_values: An entire vector can be given, fixing all parameters. Generaly you do not want to do this.
  • t=get_polyvars(rn).t: If you want to fix the time value.

kwargs

  • How parameters generally are fixed. Write e.g. n=4 to fix the parameter value of n.
DiffEqBiological.ismassactionMethod
ismassaction(network, rxidx)

Given an AbstractReactionNetwork and a reaction index, rxidx, return a boolean indicating whether the given reaction is of mass action form. For example, the reaction

2*k, 2X + 3Y --> 5Z + W

would return true, while reactions with state-dependent rates like

k*X, X + Y --> Z

would return false.

Non-allocating, returns underlying field within the reaction_network.

DiffEqBiological.jacfunMethod
jacfun(network)

Given an AbstractReactionNetwork, return a function, jac!(J,u,p,t), that evaluates the current Jacobian matrix, J, of the ODE model, $du/dt = f(u,t)$. The Jacobian matrix has entries

$J_{i,j} = \partial f_i(u,t) / \partial u_j$.

Note, for a network generated with the @min_reaction_network macro addodes! must be called first.

DiffEqBiological.jacobianexprsMethod
jacobianexprs(network)

Given an AbstractReactionNetwork, return a matrix with the ODE Jacobian expressions.

Note, for a network generated with the @min_reaction_network macro addodes! must be called first.

DiffEqBiological.jumpexprsMethod
jumpexprs(network)

Given an AbstractReactionNetwork, return a tuple of the jump rates and affects expressions.

Note, for a network generated with the @min_reaction_network macro addjumps! must be called first.

DiffEqBiological.jumpsMethod
jumps(network)

Given an AbstractReactionNetwork, return a tuple of AbstractJumps encoding a stochastical chemical kinetics representation for the reaction network.

Note, for a network generated with the @min_reaction_network macro addjumps! must be called first.

DiffEqBiological.make_hc_template!Method
make_hc_template(reaction_network)

Makes a template for solving for the given reaction networks steady states using homotopy continuation. Generally not needed as this is genrated automatically when netowkr is first solved *if one does not already exists). However, it can be sued to replace the current template, for whatever reason.

args

  • reaction_network: a reaction network.
DiffEqBiological.netstoichMethod
netstoich(network, rxidx)

Given an AbstractReactionNetwork and a reaction index, rxidx, return a vector of pairs, mapping ids of species that change numbers due to the reaction to the net stoichiometric coefficient of the species (i.e. net change in the species due to the reaction).

Allocates a new vector to store the pairs.

DiffEqBiological.noiseexprsMethod
noiseexprs(network)

Given an AbstractReactionNetwork, return a vector of the SDE noise expressions for each reaction.

Note, for a network generated with the @min_reaction_network macro addsdes! must be called first.

DiffEqBiological.noisefunMethod
noisefun(network)

Given an AbstractReactionNetwork, return a function, g(η,u,p,t), that evaluates the current noise coefficients for each reaction in the Chemical Langevin Equation representation within η.

Note, for a network generated with the @min_reaction_network macro addsdes! must be called first.

DiffEqBiological.odeexprsMethod
odeexprs(network)

Given an AbstractReactionNetwork, return a vector of the ODE expressions.

Note, for a network generated with the @min_reaction_network macro addodes! must be called first.

DiffEqBiological.odefunMethod
odefun(network)

Given an AbstractReactionNetwork, return a DiffEqBase.ODEFunction encoding an ODE model for the reaction network.

Note, for a network generated with the @min_reaction_network macro addodes! must be called first.

DiffEqBiological.oderatelawexprMethod
oderatelawexpr(network, rxidx)

Given an AbstractReactionNetwork, return the reaction rate law expression used in generated ODEs for the reaction with index rxidx. Note, for a reaction defined by

k*X*Y, X+Z --> 2X + Y

the expression that is returned will be :(k*X^2*Y*Z). For a reaction of the form

k, 2X+3Y --> Z

the expression that is returned will be :(k * (X^2/2) * (Y^3/6)).

DiffEqBiological.oderhsfunMethod
oderhsfun(network)

Given an AbstractReactionNetwork, return a function, f!(du,u,p,t), that evaluates the current value of the ODE model derivative functions, $du/dt = f(u,t)$, within du.

Note, for a network generated with the @min_reaction_network macro addodes! must be called first.

DiffEqBiological.paramjacfunMethod
paramjacfun(network)

Given an AbstractReactionNetwork, return a function, pjac(pJ,u,p,t), that evaluates the current parameter Jacobian matrix, pJ, of the ODE model, $du/dt = f(u,t)$. The parameter Jacobian matrix has entries

$pJ_{i,j} = \partial f_i(u,t) / \partial p_j$.

Note, for a network generated with the @min_reaction_network macro addodes! must be called first.

DiffEqBiological.paramsMethod
params(network)

Given an AbstractReactionNetwork, return a vector of parameter symbols.

DiffEqBiological.paramsmapMethod
paramsmap(network)

Given an AbstractReactionNetwork, return a Dictionary mapping from parameter symbol to parameter index.

DiffEqBiological.productsMethod
products(network, rxidx)

Given an AbstractReactionNetwork and a reaction index, rxidx, return a vector of symbols of species that correspond to products in the reaction. i.e. for

k*W, X + 3Y --> X + W

the returned vector would be [:X,:W].

Allocates a new vector to store the symbols.

DiffEqBiological.productstoichMethod
productstoich(network, rxidx)

Given an AbstractReactionNetwork and a reaction index, rxidx, return a vector of pairs, mapping ids of species that are products in the reaction to the corresponding stoichiometric coefficient as a product.

Allocates a new vector to store the pairs.

DiffEqBiological.productsymstoichMethod
productsymstoich(network, rxidx)

Given an AbstractReactionNetwork and a reaction index, rxidx, return a Vector of ReactantStructs, mapping the symbols of species that are products in the reaction to the corresponding stoichiometric coefficient as a product.

Non-allocating, returns underlying field within the reaction_network.

DiffEqBiological.rateexprMethod
rateexpr(network, rxidx)

Given an AbstractReactionNetwork, return the reaction rate expression for the reaction with index rxidx. Note, for a reaction defined by

k*X*Y, X+Z --> 2X + Y

the expression that is returned will be :(k*X*Y), while the rate law used in ODEs and SDEs would be k*X^2*Y*Z.

DiffEqBiological.regularjumpsMethod
regularjumps(network)

Given an AbstractReactionNetwork, return a RegularJump encoding a stochastical chemical kinetics representation of the reaction network for use in $\tau$-leaping approximations.

Note, for a network generated with the @min_reaction_network macro addjumps! must be called first.

DiffEqBiological.rxtorx_depgraphFunction
rxtorx_depgraph(network)

Given an AbstractReactionNetwork, returns a Vector{Vector{Int}} mapping each reaction index to the indices of reactions that depend on it.

DiffEqBiological.rxtospecies_depgraphMethod
rxtospecies_depgraph(network)

Given an AbstractReactionNetwork, returns a Vector{Vector{Int}} mapping each reaction index to the indices of species that depend on it.

DiffEqBiological.sdefunMethod
sdefun(network)

Given an AbstractReactionNetwork, return a DiffEqBase.SDEFunction encoding a Chemical Langevin Equation SDE model for the reaction network.

Note, for a network generated with the @min_reaction_network macro addsdes! must be called first.

DiffEqBiological.speciesMethod
species(network)

Given an AbstractReactionNetwork, return a vector of species symbols.

DiffEqBiological.speciesmapMethod
speciesmap(network)

Given an AbstractReactionNetwork, return a Dictionary mapping from species symbol to species index.

DiffEqBiological.speciestorx_depgraphMethod
speciestorx_depgraph(network)

Given an AbstractReactionNetwork, returns a Vector{Vector{Int}} mapping each species index to the indices of reactions that depend on it.

DiffEqBiological.ssaratelawexprMethod
ssaratelawexpr(network, rxidx)

Given an AbstractReactionNetwork, return the reaction rate law expression used in generated stochastic chemical kinetic model SSAs for the reaction with index rxidx. Note, for a reaction defined by

k*X*Y, X+Z --> 2X + Y

the expression that is returned will be :(k*X^2*Y*Z). For a reaction of the form

k, 2X+3Y --> Z

the expression that is returned will be :(k * binomial(X,2) * binomial(Y,3)).

DiffEqBiological.stabilityFunction
stability(solution,reaction_network, parameter_values,  t=0.)

Gives the steady stability of the steady state of a certain network. Returns true if the point is stable, false if it is not.

args

  • solution: the values of the reactants in the steady state.
  • parameter_values: the parameter values in for which the steady state was calculated.
  • reaction_network: a reaction network.
  • t (optional): the value of t, for which the stability is calculated.
DiffEqBiological.stabilityFunction
stability(solutions,reaction_network, parameter_values,  t=0.)

Gives the steady stability of the steady state of a certain network.

args

  • solutions: An array of steady states (which each in turn, is an array of reactant concentrations in the steady state).
  • parameter_values: the parameter values in for which the steady state was calculated.
  • reaction_network: a reaction network.
  • t (optional): the value of t, for which the stability is calculated.
DiffEqBiological.steady_statesMethod
steady_states(reaction_network, [solver], [parameter_values]; kwargs...)

Finds the steady states of a given reaction network, at given parameter values.

args

  • reaction_network: a reaction network.
  • solver (optional): a subtype of AbstractSteadyStateSolver that specifies how to solve the bifurcation diagram. Default is HCSteadyStateSolver (currently only solver avaiable).
  • parameter_values: a vector which specifies the parameter values for which the steady states are to be found. Not required in case the reaction network lacks steady states.

kwargs

  • potential arguments for steady state solvers (however, only current solver does not have any such arguments).
DiffEqBiological.substratesMethod
substrates(network, rxidx)

Given an AbstractReactionNetwork and a reaction index, rxidx, return a vector of symbols of species that correspond to substrates in the reaction. i.e. for

k*W, X + 3Y --> X + W

the returned vector would be [:X,:Y].

Allocates a new vector to store the symbols.

DiffEqBiological.substratestoichMethod
substratestoich(network, rxidx)

Given an AbstractReactionNetwork and a reaction index, rxidx, return a vector of pairs, mapping ids of species that serve as substrates in the reaction to the corresponding stoichiometric coefficient as a substrate.

Allocates a new vector to store the pairs.

DiffEqBiological.substratesymstoichMethod
substratesymstoich(network, rxidx)

Given an AbstractReactionNetwork and a reaction index, rxidx, return a Vector of ReactantStructs, mapping the symbols of species that serve as substrates in the reaction to the corresponding stoichiometric coefficient as a substrate.

Non-allocating, returns underlying field within the reaction_network.

DiffEqBiological.@add_constraintMacro
@add_contraint(reaction_network, constraints)

Some networks (like X ↔ Y) have an infinite number of potential steady states, depending on initial conditions. For these some additional contraint(s) is required. Typically in the form of fixed concentrations.

args

  • reaction_network: a reaction network.
  • constraints: Constraints, typically given as X + Y = C. Both parameters, varriables, and numbers can be given. However, parameters much be declared previously in the reaction network.
DiffEqBiological.@add_constraintsMacro
@add_contraints(reaction_network, constraints)

Same as @addconstraint, but a more natural way of adding several contraints such as: @addcontraints rn begin X1 + Y1 = C1 X2 + Y2 = C2 end

args

  • reaction_network: a reaction network.
  • constraints: Constraints, typically given as X + Y = C. Both parameters, varriables, and numbers can be given. However, parameters much be declared previously in the reaction network.
DiffEqBiological.@add_hc_templateMacro
@make_hc_template(reaction_network)

Calls the "add_hc_template function on the given reaction network."

args

  • reaction_network: a reaction network.
DiffEqBiological.@make_hc_templateMacro
@make_hc_template(reaction_network)

Calls the "make_hc_template function on the given reaction network."

args

  • reaction_network: a reaction network.
Base.:==Method
==(rn1::DiffEqBase.AbstractReactionNetwork, rn2::DiffEqBase.AbstractReactionNetwork)

Tests whether the underlying species symbols, parameter symbols and reactions are the same in the two networks. Ignores order network components were defined, so the integer id of any individual species/parameters/reactions may be different between the two networks. Does not currently account for different reaction definitions, so "k, X+Y –> Y + Z" will not be the same as "k, Y+X –> Y + Z"

DiffEqBiological.clean_subtractionsMethod

clean_subtractions(ex::Expr) Replace additions of negative terms with subtractions. This is a fairly stupid function which is designed for a specific problem with reaction networks. It is neither recursive nor very general. Return :: cleaned out expression

From Latexify.jl with permission: see

DiffEqBiological.netsymstoichMethod
netsymstoich(network, rxidx)

Given an AbstractReactionNetwork and a reaction index, rxidx, return a Vector of ReactantStructs, mapping the symbols of species that change numbers due to the reaction to the net stoichiometric coefficient of the species (i.e. net change in the species due to the reaction).

Non-allocating, returns underlying field within the reaction_network.