Macros (pysb.macros)

A collection of generally useful modeling macros.

These macros are written to be as generic and reusable as possible, serving as a collection of best practices and implementation ideas. They conform to the following general guidelines:

  • All components created by the macro are implicitly added to the current model and explicitly returned in a ComponentSet.

  • Parameters may be passed as Parameter or Expression objects, or as plain numbers for which Parameter objects will be automatically created using an appropriate naming convention.

  • Arguments which accept a MonomerPattern should also accept Monomers, which are to be interpreted as MonomerPatterns on that Monomer with an empty condition list. This is typically implemented by having the macro apply the “call” (parentheses) operator to the argument with an empty argument list and using the resulting value instead of the original argument when creating Rules, e.g. arg = arg(). Calling a Monomer will return a MonomerPattern, and calling a MonomerPattern will return a copy of itself, so calling either is guaranteed to return a MonomerPattern.

The _macro_rule helper function contains much of the logic needed to follow these guidelines. Every macro in this module either uses _macro_rule directly or calls another macro which does.

Another useful function is _verify_sites which will raise an exception if a Monomer or MonomerPattern does not possess every one of a given list of sites. This can be used to trigger such errors up front rather than letting an exception occur at the point where the macro tries to use the invalid site in a pattern, which can be harder for the caller to debug.

pysb.macros.assemble_chain_sequential_base(base, basesite, subunit, site1, site2, max_size, ktable, comp=1)[source]

Generate rules to assemble a homomeric chain sequentially onto a base complex (only the subunit creates repeating chain, not the base).

The chain species are created by sequential addition of subunit monomers. The chain structure is defined by the pore_species_base macro.

Parameters:
baseMonomer or MonomerPattern

The base complex to which the chain is attached.

basesitestring

The name of the site on the complex to which chain attaches.

subunitMonomer or MonomerPattern

The subunit of which the chain is composed.

site1, site2string

The names of the sites where one copy of subunit binds to the next; the first will also be the site where the first subunit binds the base.

max_sizeinteger

The maximum number of subunits in the chain.

ktablelist of lists of Parameters or numbers

Table of forward and reverse rate constants for the assembly steps. The outer list must be of length max_size + 1, and the inner lists must all be of length 2. In the outer list, the first element corresponds to the first assembly step in which the complex binds the first subunit. The next corresponds to a bound subunit binding to form a 2-subunit complex, and the last element corresponds to the final step in which the max_size`th subunit is added. Each inner list contains the forward and reverse rate constants (in that order) for the corresponding assembly reaction, and each of these pairs must comprise solely Parameter objects or solely numbers (never one of each). If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on `subunit, site1, site2 and the chain sizes and these parameters will be included at the end of the returned component list.

compoptional; a ComplexPattern to which the base molecule is attached.

Examples

Assemble a three-membered chain by sequential addition of monomers to a base, which is in turn attached to a complex, with the same forward/reverse rates for monomer-monomer and monomer-dimer interactions:

Model()
Monomer('Base', ['b1', 'b2'])
Monomer('Unit', ['p1', 'p2'])
Monomer('Complex1', ['s1'])
Monomer('Complex2', ['s1', s2'])
assemble_chain_sequential(Base(b2=ANY), 'b1', Unit, 'p1', 'p2', 3, [[1e-4, 1e-1]] * 2, Complex1(s1=ANY) % Complex2(s1=ANY, s2=ANY))

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('Base', ['b1', 'b2'])
Monomer('Base', ['b1', 'b2'])
>>> Monomer('Unit', ['p1', 'p2'])
Monomer('Unit', ['p1', 'p2'])
>>> Monomer('Complex1', ['s1'])
Monomer('Complex1', ['s1'])
>>> Monomer('Complex2', ['s1', 's2'])
Monomer('Complex2', ['s1', 's2'])
>>> assemble_chain_sequential_base(Base(b2=ANY), 'b1', Unit, 'p1', 'p2', 3, [[1e-4, 1e-1]] * 2, Complex1(s1=ANY) % Complex2(s1=ANY, s2=ANY)) 
ComponentSet([
 Rule('assemble_chain_sequential_base_Unit_2', Unit(p1=None, p2=None) + Complex1(s1=ANY) % Complex2(s1=ANY, s2=ANY) % Base(b1=1, b2=ANY) % Unit(p1=1, p2=None) | Complex1(s1=ANY) % Complex2(s1=ANY, s2=ANY) % Base(b1=1, b2=ANY) % Unit(p1=1, p2=2) % Unit(p1=2, p2=None), assemble_chain_sequential_base_Unit_2_kf, assemble_chain_sequential_base_Unit_2_kr),
 Parameter('assemble_chain_sequential_base_Unit_2_kf', 0.0001),
 Parameter('assemble_chain_sequential_base_Unit_2_kr', 0.1),
 Rule('assemble_chain_sequential_base_Unit_3', Unit(p1=None, p2=None) + Complex1(s1=ANY) % Complex2(s1=ANY, s2=ANY) % Base(b1=1, b2=ANY) % Unit(p1=1, p2=2) % Unit(p1=2, p2=None) | MatchOnce(Complex1(s1=ANY) % Complex2(s1=ANY, s2=ANY) % Base(b1=1, b2=ANY) % Unit(p1=1, p2=2) % Unit(p1=2, p2=3) % Unit(p1=3, p2=None)), assemble_chain_sequential_base_Unit_3_kf, assemble_chain_sequential_base_Unit_3_kr),
 Parameter('assemble_chain_sequential_base_Unit_3_kf', 0.0001),
 Parameter('assemble_chain_sequential_base_Unit_3_kr', 0.1),
 ])
pysb.macros.assemble_pore_sequential(subunit, site1, site2, max_size, ktable)[source]

Generate rules to assemble a circular homomeric pore sequentially.

The pore species are created by sequential addition of subunit monomers, i.e. larger oligomeric species never fuse together. The pore structure is defined by the pore_species macro.

Parameters:
subunitMonomer or MonomerPattern

The subunit of which the pore is composed.

site1, site2string

The names of the sites where one copy of subunit binds to the next.

max_sizeinteger

The maximum number of subunits in the pore.

ktablelist of lists of Parameters or numbers

Table of forward and reverse rate constants for the assembly steps. The outer list must be of length max_size - 1, and the inner lists must all be of length 2. In the outer list, the first element corresponds to the first assembly step in which two monomeric subunits bind to form a 2-subunit complex, and the last element corresponds to the final step in which the max_size`th subunit is added. Each inner list contains the forward and reverse rate constants (in that order) for the corresponding assembly reaction, and each of these pairs must comprise solely Parameter objects or solely numbers (never one of each). If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on `subunit, site1, site2 and the pore sizes and these parameters will be included at the end of the returned component list.

Examples

Assemble a three-membered pore by sequential addition of monomers, with the same forward/reverse rates for monomer-monomer and monomer-dimer interactions:

Model()
Monomer('Unit', ['p1', 'p2'])
assemble_pore_sequential(Unit, 'p1', 'p2', 3, [[1e-4, 1e-1]] * 2)

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('Unit', ['p1', 'p2'])
Monomer('Unit', ['p1', 'p2'])
>>> assemble_pore_sequential(Unit, 'p1', 'p2', 3, [[1e-4, 1e-1]] * 2) 
ComponentSet([
 Rule('assemble_pore_sequential_Unit_2',
      Unit(p1=None, p2=None) + Unit(p1=None, p2=None) |
          Unit(p1=None, p2=1) % Unit(p1=1, p2=None),
      assemble_pore_sequential_Unit_2_kf,
      assemble_pore_sequential_Unit_2_kr),
 Parameter('assemble_pore_sequential_Unit_2_kf', 0.0001),
 Parameter('assemble_pore_sequential_Unit_2_kr', 0.1),
 Rule('assemble_pore_sequential_Unit_3',
      Unit(p1=None, p2=None) + Unit(p1=None, p2=1) % Unit(p1=1, p2=None) |
          MatchOnce(Unit(p1=3, p2=1) % Unit(p1=1, p2=2) % Unit(p1=2, p2=3)),
      assemble_pore_sequential_Unit_3_kf,
      assemble_pore_sequential_Unit_3_kr),
 Parameter('assemble_pore_sequential_Unit_3_kf', 0.0001),
 Parameter('assemble_pore_sequential_Unit_3_kr', 0.1),
 ])
pysb.macros.bind(s1, site1, s2, site2, klist)[source]

Generate the reversible binding reaction S1 + S2 | S1:S2.

Parameters:
s1, s2Monomer or MonomerPattern

Monomers participating in the binding reaction.

site1, site2string

The names of the sites on s1 and s2 used for binding.

klistlist of 2 Parameters or list of 2 numbers

Forward and reverse rate constants (in that order). If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on the names and states of S1 and S2 and these parameters will be included at the end of the returned component list.

Returns:
componentsComponentSet

The generated components. Contains the bidirectional binding Rule and optionally two Parameters if klist was given as numbers.

Examples

Binding between A and B:

Model()
Monomer('A', ['x'])
Monomer('B', ['y'])
bind(A, 'x', B, 'y', [1e-4, 1e-1])

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('A', ['x'])
Monomer('A', ['x'])
>>> Monomer('B', ['y'])
Monomer('B', ['y'])
>>> bind(A, 'x', B, 'y', [1e-4, 1e-1]) 
ComponentSet([
 Rule('bind_A_B', A(x=None) + B(y=None) | A(x=1) % B(y=1), bind_A_B_kf, bind_A_B_kr),
 Parameter('bind_A_B_kf', 0.0001),
 Parameter('bind_A_B_kr', 0.1),
 ])
pysb.macros.bind_complex(s1, site1, s2, site2, klist, m1=None, m2=None)[source]

Generate the reversible binding reaction S1 + S2 | S1:S2, with optional complexes attached to either S1 (C1:S1 + S2 | C1:S1:S2), S2 (S1 + C2:S2 | C2:S2:S1), or both (C1:S1 + C2:S2 | C1:S1:S2:C2).

Parameters:
s1, s2Monomer, MonomerPattern, or ComplexPattern

Monomers or complexes participating in the binding reaction.

site1, site2string

The names of the sites on s1 and s2 used for binding.

klistlist of 2 Parameters or list of 2 numbers

Forward and reverse rate constants (in that order). If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on the names and states of S1 and S2 and these parameters will be included at the end of the returned component list.

m1, m2Monomer or MonomerPattern

If s1 or s2 binding site is present in multiple monomers within a complex, the specific monomer desired for binding must be specified.

Returns:
componentsComponentSet

The generated components. Contains the bidirectional binding Rule and optionally two Parameters if klist was given as numbers.

Examples

Binding between A:B and C:D:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('A', ['a', 'b'])
Monomer('A', ['a', 'b'])
>>> Monomer('B', ['c', 'd'])
Monomer('B', ['c', 'd'])
>>> Monomer('C', ['e', 'f'])
Monomer('C', ['e', 'f'])
>>> Monomer('D', ['g', 'h'])
Monomer('D', ['g', 'h'])
>>> bind_complex(A(a=1) % B(c=1), 'b', C(e=2) % D(g=2), 'h', [1e-4,             1e-1]) 
ComponentSet([
Rule('bind_AB_DC', A(a=1, b=None) % B(c=1) + D(g=3, h=None) % C(e=3)
  | A(a=1, b=50) % B(c=1) % D(g=3, h=50) % C(e=3), bind_AB_DC_kf,
  bind_AB_DC_kr),
Parameter('bind_AB_DC_kf', 0.0001),
Parameter('bind_AB_DC_kr', 0.1),
])

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('A', ['a', 'b'])
Monomer('A', ['a', 'b'])
>>> Monomer('B', ['c', 'd'])
Monomer('B', ['c', 'd'])
>>> Monomer('C', ['e', 'f'])
Monomer('C', ['e', 'f'])
>>> Monomer('D', ['g', 'h'])
Monomer('D', ['g', 'h'])
>>> bind(A, 'a', B, 'c', [1e4, 1e-1]) 
ComponentSet([
Rule('bind_A_B',
  A(a=None) + B(c=None) | A(a=1) % B(c=1),
  bind_A_B_kf, bind_A_B_kr),
Parameter('bind_A_B_kf', 10000.0),
Parameter('bind_A_B_kr', 0.1),
])
>>> bind(C, 'e', D, 'g', [1e4, 1e-1]) 
ComponentSet([
Rule('bind_C_D',
  C(e=None) + D(g=None) | C(e=1) % D(g=1),
  bind_C_D_kf, bind_C_D_kr),
Parameter('bind_C_D_kf', 10000.0),
Parameter('bind_C_D_kr', 0.1),
])
>>> bind_complex(A(a=1) % B(c=1), 'b', C(e=2) % D(g=2), 'h', [1e-4,             1e-1]) 
ComponentSet([
Rule('bind_AB_DC',
  A(a=1, b=None) % B(c=1) + D(g=3, h=None) % C(e=3) | A(a=1,
  b=50) % B(c=1) % D(g=3, h=50) % C(e=3),
  bind_AB_DC_kf, bind_AB_DC_kr),
Parameter('bind_AB_DC_kf', 0.0001),
Parameter('bind_AB_DC_kr', 0.1),
])
pysb.macros.bind_table(bindtable, row_site, col_site, kf=None)[source]

Generate a table of reversible binding reactions.

Given two lists of species R and C, calls the bind macro on each pairwise combination (R[i], C[j]). The species lists and the parameter values are passed as a list of lists (i.e. a table) with elements of R passed as the “row headers”, elements of C as the “column headers”, and forward / reverse rate pairs (in that order) as tuples in the “cells”. For example with two elements in each of R and C, the table would appear as follows (note that the first row has one fewer element than the subsequent rows):

[[              C1,           C2],
 [R1, (1e-4, 1e-1), (2e-4, 2e-1)],
 [R2, (3e-4, 3e-1), (4e-4, 4e-1)]]

Each parameter tuple may contain Parameters or numbers. If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on the names and states of the relevant species and these parameters will be included at the end of the returned component list. To omit any individual reaction, pass None in place of the corresponding parameter tuple.

Alternately, single kd values (dissociation constant, kr/kf) may be specified instead of (kf, kr) tuples. If kds are used, a single shared kf Parameter or number must be passed as an extra kf argument. kr values for each binding reaction will be calculated as kd*kf. It is important to remember that the forward rate constant is a single parameter shared across the entire bind table, as this may have implications for parameter fitting.

Parameters:
bindtablelist of lists

Table of reactants and rates, as described above.

row_site, col_sitestring

The names of the sites on the elements of R and C, respectively, used for binding.

kfParameter or number, optional

If the “cells” in bindtable are given as single kd values, this is the shared kf used to calculate the kr values.

Returns:
componentsComponentSet

The generated components. Contains the bidirectional binding Rules and optionally the Parameters for any parameters given as numbers.

Examples

Binding table for two species types (R and C), each with two members:

Model()
Monomer('R1', ['x'])
Monomer('R2', ['x'])
Monomer('C1', ['y'])
Monomer('C2', ['y'])
bind_table([[               C1,           C2],
            [R1,  (1e-4, 1e-1),  (2e-4, 2e-1)],
            [R2,  (3e-4, 3e-1),         None]],
           'x', 'y')

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('R1', ['x'])
Monomer('R1', ['x'])
>>> Monomer('R2', ['x'])
Monomer('R2', ['x'])
>>> Monomer('C1', ['y'])
Monomer('C1', ['y'])
>>> Monomer('C2', ['y'])
Monomer('C2', ['y'])
>>> bind_table([[               C1,           C2],
...             [R1,  (1e-4, 1e-1),  (2e-4, 2e-1)],
...             [R2,  (3e-4, 3e-1),         None]],
...            'x', 'y') 
ComponentSet([
 Rule('bind_R1_C1', R1(x=None) + C1(y=None) | R1(x=1) % C1(y=1),
     bind_R1_C1_kf, bind_R1_C1_kr),
 Parameter('bind_R1_C1_kf', 0.0001),
 Parameter('bind_R1_C1_kr', 0.1),
 Rule('bind_R1_C2', R1(x=None) + C2(y=None) | R1(x=1) % C2(y=1),
     bind_R1_C2_kf, bind_R1_C2_kr),
 Parameter('bind_R1_C2_kf', 0.0002),
 Parameter('bind_R1_C2_kr', 0.2),
 Rule('bind_R2_C1', R2(x=None) + C1(y=None) | R2(x=1) % C1(y=1),
     bind_R2_C1_kf, bind_R2_C1_kr),
 Parameter('bind_R2_C1_kf', 0.0003),
 Parameter('bind_R2_C1_kr', 0.3),
 ])
pysb.macros.bind_table_complex(bindtable, row_site, col_site, m1=None, m2=None, kf=None)[source]

Generate a table of reversible binding reactions when either the row or column species (or both) have a complex bound to them.

Given two lists of species R and C (which can be complexes or monomers), calls the bind_complex macro on each pairwise combination (R[i], C[j]). The species lists and the parameter values are passed as a list of lists (i.e. a table) with elements of R passed as the “row headers”, elements of C as the “column headers”, and forward / reverse rate pairs (in that order) as tuples in the “cells”. For example with two elements in each of R and C, the table would appear as follows (note that the first row has one fewer element than the subsequent rows):

[[              C1,           C2],
 [R1, (1e-4, 1e-1), (2e-4, 2e-1)],
 [R2, (3e-4, 3e-1), (4e-4, 4e-1)]]

Each parameter tuple may contain Parameters or numbers. If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on the names and states of the relevant species and these parameters will be included at the end of the returned component list. To omit any individual reaction, pass None in place of the corresponding parameter tuple.

Alternately, single kd values (dissociation constant, kr/kf) may be specified instead of (kf, kr) tuples. If kds are used, a single shared kf Parameter or number must be passed as an extra kf argument. kr values for each binding reaction will be calculated as kd*kf. It is important to remember that the forward rate constant is a single parameter shared across the entire bind table, as this may have implications for parameter fitting.

Parameters:
bindtablelist of lists

Table of reactants and rates, as described above.

row_site, col_sitestring

The names of the sites on the elements of R and C, respectively, used for binding.

m1Monomer or MonomerPattern, optional

Monomer in row complex for binding. Must be specified if there are multiple monomers that have the row_site within a complex.

m2Monomer or MonomerPattern, optional

Monomer in column complex for binding. Must be specified if there are multiple monomers that have the col_site within a complex.

kfParameter or number, optional

If the “cells” in bindtable are given as single kd values, this is the shared kf used to calculate the kr values.

Returns:
componentsComponentSet

The generated components. Contains the bidirectional binding Rules and optionally the Parameters for any parameters given as numbers.

Examples

Binding table for two species types (R and C, which can be complexes or monomers):

Model()
Monomer('R1', ['x', 'c1'])
Monomer('R2', ['x', 'c1'])
Monomer('C1', ['y', 'c2'])
Monomer('C2', ['y', 'c2'])
bind(C1(y=None), 'c2', C1(y=None), 'c2', (1e-3, 1e-2))
bind(R1(x=None), 'c1', R2(x=None), 'c1', (1e-3, 1e-2))
bind_table_complex([[              C1(c2=1, y=None)%C1(c2=1),            C2],
                   [R1()%R2(),  (1e-4, 1e-1),  (2e-4, 2e-1)],
                   [R2,     (3e-4, 3e-1),         None]],
                   'x', 'y', m1=R1(), m2=C1(y=None, c2=1))

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('R1', ['x', 'c1'])
Monomer('R1', ['x', 'c1'])
>>> Monomer('R2', ['x', 'c1'])
Monomer('R2', ['x', 'c1'])
>>> Monomer('C1', ['y', 'c2'])
Monomer('C1', ['y', 'c2'])
>>> Monomer('C2', ['y', 'c2'])
Monomer('C2', ['y', 'c2'])
>>> bind(C1(y=None), 'c2', C1(y=None), 'c2', (1e-3, 1e-2)) 
ComponentSet([
 Rule('bind_C1_C1', C1(y=None, c2=None) + C1(y=None, c2=None) | C1(y=None, c2=1) % C1(y=None, c2=1), bind_C1_C1_kf, bind_C1_C1_kr),
 Parameter('bind_C1_C1_kf', 0.001),
 Parameter('bind_C1_C1_kr', 0.01),
 ])
>>> bind(R1(x=None), 'c1', R2(x=None), 'c1', (1e-3, 1e-2)) 
ComponentSet([
 Rule('bind_R1_R2', R1(x=None, c1=None) + R2(x=None, c1=None) | R1(x=None, c1=1) % R2(x=None, c1=1), bind_R1_R2_kf, bind_R1_R2_kr),
 Parameter('bind_R1_R2_kf', 0.001),
 Parameter('bind_R1_R2_kr', 0.01),
 ])
>>> bind_table_complex([[               C1(c2=1, y=None)%C1(c2=1),           C2],
...                      [R1()%R2(),      (1e-4, 1e-1),                        (2e-4, 2e-1)],
...                       [R2,             (3e-4, 3e-1),                        None]],
...                       'x', 'y', m1=R1(), m2=C1(y=None, c2=1)) 
ComponentSet([
Rule('bind_R1R2_C1C1', R1(x=None) % R2() + C1(y=None, c2=1) % C1(c2=1) | R1(x=50) % R2() % C1(y=50, c2=1) % C1(c2=1), bind_R1R2_C1C1_kf, bind_R1R2_C1C1_kr),
Parameter('bind_R1R2_C1C1_kf', 0.0001),
Parameter('bind_R1R2_C1C1_kr', 0.1),
Rule('bind_R1R2_C2', R1(x=None) % R2() + C2(y=None) | R1(x=50) % R2() % C2(y=50), bind_R1R2_C2_kf, bind_R1R2_C2_kr),
Parameter('bind_R1R2_C2_kf', 0.0002),
Parameter('bind_R1R2_C2_kr', 0.2),
Rule('bind_C1C1_R2', C1(y=None, c2=1) % C1(c2=1) + R2(x=None) | C1(y=50, c2=1) % C1(c2=1) % R2(x=50), bind_C1C1_R2_kf, bind_C1C1_R2_kr),
Parameter('bind_C1C1_R2_kf', 0.0003),
Parameter('bind_C1C1_R2_kr', 0.3),
 ])
pysb.macros.catalyze(enzyme, e_site, substrate, s_site, product, klist)[source]

Generate the two-step catalytic reaction E + S | E:S >> E + P.

Parameters:
enzyme, substrate, productMonomer or MonomerPattern

E, S and P in the above reaction.

e_site, s_sitestring

The names of the sites on enzyme and substrate (respectively) where they bind each other to form the E:S complex.

klistlist of 3 Parameters or list of 3 numbers

Forward, reverse and catalytic rate constants (in that order). If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on the names and states of enzyme, substrate and product and these parameters will be included at the end of the returned component list.

Returns:
componentsComponentSet

The generated components. Contains two Rules (bidirectional complex formation and unidirectional product dissociation), and optionally three Parameters if klist was given as plain numbers.

Notes

When passing a MonomerPattern for enzyme or substrate, do not include e_site or s_site in the respective patterns. The macro will handle this.

Examples

Using distinct Monomers for substrate and product:

Model()
Monomer('E', ['b'])
Monomer('S', ['b'])
Monomer('P')
catalyze(E(), 'b', S(), 'b', P(), (1e-4, 1e-1, 1))

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('E', ['b'])
Monomer('E', ['b'])
>>> Monomer('S', ['b'])
Monomer('S', ['b'])
>>> Monomer('P')
Monomer('P')
>>> catalyze(E(), 'b', S(), 'b', P(), (1e-4, 1e-1, 1)) 
ComponentSet([
 Rule('bind_E_S_to_ES', E(b=None) + S(b=None) | E(b=1) % S(b=1),
     bind_E_S_to_ES_kf, bind_E_S_to_ES_kr),
 Parameter('bind_E_S_to_ES_kf', 0.0001),
 Parameter('bind_E_S_to_ES_kr', 0.1),
 Rule('catalyze_ES_to_E_P', E(b=1) % S(b=1) >> E(b=None) + P(),
     catalyze_ES_to_E_P_kc),
 Parameter('catalyze_ES_to_E_P_kc', 1.0),
 ])

Using a single Monomer for substrate and product with a state change:

Monomer('Kinase', ['b'])
Monomer('Substrate', ['b', 'y'], {'y': ('U', 'P')})
catalyze(Kinase(), 'b', Substrate(y='U'), 'b', Substrate(y='P'),
         (1e-4, 1e-1, 1))

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('Kinase', ['b'])
Monomer('Kinase', ['b'])
>>> Monomer('Substrate', ['b', 'y'], {'y': ('U', 'P')})
Monomer('Substrate', ['b', 'y'], {'y': ('U', 'P')})
>>> catalyze(Kinase(), 'b', Substrate(y='U'), 'b', Substrate(y='P'), (1e-4, 1e-1, 1)) 
ComponentSet([
 Rule('bind_Kinase_SubstrateU_to_KinaseSubstrateU',
     Kinase(b=None) + Substrate(b=None, y='U') | Kinase(b=1) % Substrate(b=1, y='U'),
     bind_Kinase_SubstrateU_to_KinaseSubstrateU_kf,
     bind_Kinase_SubstrateU_to_KinaseSubstrateU_kr),
 Parameter('bind_Kinase_SubstrateU_to_KinaseSubstrateU_kf', 0.0001),
 Parameter('bind_Kinase_SubstrateU_to_KinaseSubstrateU_kr', 0.1),
 Rule('catalyze_KinaseSubstrateU_to_Kinase_SubstrateP',
      Kinase(b=1) % Substrate(b=1, y='U') >> Kinase(b=None) + Substrate(b=None, y='P'),
      catalyze_KinaseSubstrateU_to_Kinase_SubstrateP_kc),
 Parameter('catalyze_KinaseSubstrateU_to_Kinase_SubstrateP_kc', 1.0),
 ])
pysb.macros.catalyze_complex(enzyme, e_site, substrate, s_site, product, klist, m1=None, m2=None)[source]

Generate the two-step catalytic reaction E + S | E:S >> E + P, while allowing complexes to serve as enzyme, substrate and/or product.

E:S1 + S:S2 | E:S1:S:S2 >> E:S1 + P:S2

Parameters:
enzyme, substrate, productMonomer, MonomerPattern, or ComplexPattern
Monomers or complexes participating in the binding reaction.
e_site, s_sitestring
The names of the sites on `enzyme` and `substrate` (respectively) where
they bind each other to form the E:S complex.
klistlist of 3 Parameters or list of 3 numbers
Forward, reverse and catalytic rate constants (in that order). If
Parameters are passed, they will be used directly in the generated
Rules. If numbers are passed, Parameters will be created with
automatically generated names based on the names and states of enzyme,
substrate and product and these parameters will be included at the end
of the returned component list.
m1, m2Monomer or MonomerPattern
If enzyme or substrate binding site is present in multiple monomers
within a complex, the specific monomer desired for binding must be specified.
Returns:
componentsComponentSet
The generated components. Contains the bidirectional binding Rule
and optionally three Parameters if klist was given as numbers.
pysb.macros.catalyze_one_step(enzyme, substrate, product, kf)[source]

Generate the one-step catalytic reaction E + S >> E + P.

Parameters:
enzyme, substrate, productMonomer or MonomerPattern

E, S and P in the above reaction.

kfa Parameter or a number

Forward rate constant for the reaction. If a Parameter is passed, it will be used directly in the generated Rules. If a number is passed, a Parameter will be created with an automatically generated name based on the names and states of the enzyme, substrate and product and this parameter will be included at the end of the returned component list.

Returns:
componentsComponentSet

The generated components. Contains the unidirectional reaction Rule and optionally the forward rate Parameter if klist was given as a number.

Notes

In this macro, there is no direct binding between enzyme and substrate, so binding sites do not have to be specified. This represents an approximation for the case when the enzyme is operating in its linear range. However, if catalysis is nevertheless contingent on the enzyme or substrate being unbound on some site, then that information must be encoded in the MonomerPattern for the enzyme or substrate. See the examples, below.

Examples

Convert S to P by E:

Model()
Monomer('E', ['b'])
Monomer('S', ['b'])
Monomer('P')
catalyze_one_step(E, S, P, 1e-4)

If the ability of the enzyme E to catalyze this reaction is dependent on the site ‘b’ of E being unbound, then this macro must be called as

catalyze_one_step(E(b=None), S, P, 1e-4)

and similarly if the substrate or product must be unbound.

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('E', ['b'])
Monomer('E', ['b'])
>>> Monomer('S', ['b'])
Monomer('S', ['b'])
>>> Monomer('P')
Monomer('P')
>>> catalyze_one_step(E, S, P, 1e-4) 
ComponentSet([
 Rule('one_step_E_S_to_E_P', E() + S() >> E() + P(), one_step_E_S_to_E_P_kf),
 Parameter('one_step_E_S_to_E_P_kf', 0.0001),
 ])
pysb.macros.catalyze_one_step_reversible(enzyme, substrate, product, klist)[source]

Create fwd and reverse rules for catalysis of the form:

E + S -> E + P
    P -> S 
Parameters:
enzyme, substrate, productMonomer or MonomerPattern

E, S and P in the above reactions.

klistlist of 2 Parameters or list of 2 numbers

A list containing the rate constant for catalysis and the rate constant for the conversion of product back to substrate (in that order). If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on the names and states of S1 and S2 and these parameters will be included at the end of the returned component list.

Returns:
componentsComponentSet

The generated components. Contains two rules (the single-step catalysis rule and the product reversion rule) and optionally the two generated Parameter objects if klist was given as numbers.

Notes

Calls the macro catalyze_one_step to generate the catalysis rule.

Examples

One-step, pseudo-first order conversion of S to P by E:

Model()
Monomer('E', ['b'])
Monomer('S', ['b'])
Monomer('P')
catalyze_one_step_reversible(E, S, P, [1e-1, 1e-4])

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('E', ['b'])
Monomer('E', ['b'])
>>> Monomer('S', ['b'])
Monomer('S', ['b'])
>>> Monomer('P')
Monomer('P')
>>> catalyze_one_step_reversible(E, S, P, [1e-1, 1e-4]) 
ComponentSet([
 Rule('one_step_E_S_to_E_P', E() + S() >> E() + P(), one_step_E_S_to_E_P_kf),
 Parameter('one_step_E_S_to_E_P_kf', 0.1),
 Rule('reverse_P_to_S', P() >> S(), reverse_P_to_S_kr),
 Parameter('reverse_P_to_S_kr', 0.0001),
 ])
pysb.macros.catalyze_state(enzyme, e_site, substrate, s_site, mod_site, state1, state2, klist)[source]

Generate the two-step catalytic reaction E + S | E:S >> E + P. A wrapper around catalyze() with a signature specifying the state change of the substrate resulting from catalysis.

Parameters:
enzymeMonomer or MonomerPattern

E in the above reaction.

substrateMonomer or MonomerPattern

S and P in the above reaction. The product species is assumed to be identical to the substrate species in all respects except the state of the modification site. The state of the modification site should not be specified in the MonomerPattern for the substrate.

e_site, s_sitestring

The names of the sites on enzyme and substrate (respectively) where they bind each other to form the E:S complex.

mod_sitestring

The name of the site on the substrate that is modified by catalysis.

state1, state2strings

The states of the modification site (mod_site) on the substrate before (state1) and after (state2) catalysis.

klistlist of 3 Parameters or list of 3 numbers

Forward, reverse and catalytic rate constants (in that order). If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on the names and states of enzyme, substrate and product and these parameters will be included at the end of the returned component list.

Returns:
componentsComponentSet

The generated components. Contains two Rules (bidirectional complex formation and unidirectional product dissociation), and optionally three Parameters if klist was given as plain numbers.

Notes

When passing a MonomerPattern for enzyme or substrate, do not include e_site or s_site in the respective patterns. In addition, do not include the state of the modification site on the substrate. The macro will handle this.

Examples

Using a single Monomer for substrate and product with a state change:

Monomer('Kinase', ['b'])
Monomer('Substrate', ['b', 'y'], {'y': ('U', 'P')})
catalyze_state(Kinase, 'b', Substrate, 'b', 'y', 'U', 'P',
         (1e-4, 1e-1, 1))

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('Kinase', ['b'])
Monomer('Kinase', ['b'])
>>> Monomer('Substrate', ['b', 'y'], {'y': ('U', 'P')})
Monomer('Substrate', ['b', 'y'], {'y': ('U', 'P')})
>>> catalyze_state(Kinase, 'b', Substrate, 'b', 'y', 'U', 'P', (1e-4, 1e-1, 1)) 
ComponentSet([
 Rule('bind_Kinase_SubstrateU_to_KinaseSubstrateU',
     Kinase(b=None) + Substrate(b=None, y='U') | Kinase(b=1) % Substrate(b=1, y='U'),
     bind_Kinase_SubstrateU_to_KinaseSubstrateU_kf,
     bind_Kinase_SubstrateU_to_KinaseSubstrateU_kr),
 Parameter('bind_Kinase_SubstrateU_to_KinaseSubstrateU_kf', 0.0001),
 Parameter('bind_Kinase_SubstrateU_to_KinaseSubstrateU_kr', 0.1),
 Rule('catalyze_KinaseSubstrateU_to_Kinase_SubstrateP',
     Kinase(b=1) % Substrate(b=1, y='U') >> Kinase(b=None) + Substrate(b=None, y='P'),
     catalyze_KinaseSubstrateU_to_Kinase_SubstrateP_kc),
 Parameter('catalyze_KinaseSubstrateU_to_Kinase_SubstrateP_kc', 1.0),
 ])
pysb.macros.degrade(species, kdeg)[source]

Generate a reaction which degrades a species.

Note that species is not required to be “concrete”.

Parameters:
speciesMonomer, MonomerPattern or ComplexPattern

The species to synthesize. If a Monomer, sites are considered as unbound and in their default state. If a pattern, must be concrete.

kdegParameters or number

Degradation rate. If a Parameter is passed, it will be used directly in the generated Rule. If a number is passed, a Parameter will be created with an automatically generated name based on the names and site states of the components of species and this parameter will be included at the end of the returned component list.

Returns:
componentsComponentSet

The generated components. Contains the unidirectional degradation Rule and optionally a Parameter if ksynth was given as a number.

Examples

Degrade all B, even bound species:

Model()
Monomer('B', ['x'])
degrade(B(), 1e-6)

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('B', ['x'])
Monomer('B', ['x'])
>>> degrade(B(), 1e-6) 
ComponentSet([
 Rule('degrade_B', B() >> None, degrade_B_k),
 Parameter('degrade_B_k', 1e-06),
 ])
pysb.macros.drug_binding(drug, d_site, substrate, s_site, t_action, klist)[source]

Generate the reversible binding reaction DRUG + SUBSTRATE | DRUG:SUBSTRATE that only gets triggered when the simulation reaches the time point t_action. The idea of this macro is to mimic experimental settings when a reaction is started and later on some kind of perturbation is added to the system.

Warning

This macro only works when a model is simulated using a deterministic simulator.

Parameters:
drug, substrate: Monomer or MonomerPattern

Monomers participating in the binding reaction.

d_site, s_site: string

The names of the sites on s1 and s2 used for binding.

t_action: float

Time of the simulation at which the drug is added

klist: list of 2 Parameters or list of 2 numbers

Forward and reverse rate constants (in that order). If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on the names and states of S1 and S2 and these parameters will be included at the end of the returned component list.

Returns:
componentsComponentSet

The generated components. Contains the bidirectional binding Rule, the time monomer, Parameter rate of time creation, Rule to simulate passing of time, time Observable, two Expression rates that take into account when the interaction between the drug and that substrate start to occur and optionally two Parameters if klist was given as numbers as numbers

Examples

Binding between drug and substrate::

Model() Monomer(‘drug’, [‘b’]) Monomer(‘substrate’, [‘b’]) drug_binding(drug(), ‘b’, substrate(), ‘b’, 10, [2,4])

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('drug', ['b'])
Monomer('drug', ['b'])
>>> Monomer('substrate', ['b'])
Monomer('substrate', ['b'])
>>> drug_binding(drug(), 'b', substrate(), 'b', 10, [0.1, 0.01])
ComponentSet([
 Rule('bind_drug_substrate_to_drugsubstrate', drug(b=None) + substrate(b=None) | drug(b=1) % substrate(b=1), kf_expr_drug_substrate, kr_expr_drug_substrate),
 Parameter('kf_drug_substrate', 0.1),
 Parameter('kr_drug_substrate', 0.01),
 Rule('synthesize___t', None >> __t(), __k_t),
 Monomer('__t'),
 Parameter('__k_t', 1.0),
 Observable('t', __t()),
 Expression('kf_expr_drug_substrate', Piecewise((kf_drug_substrate, t > 10), (0, True))),
 Expression('kr_expr_drug_substrate', Piecewise((kr_drug_substrate, t > 10), (0, True))),
 ])
pysb.macros.equilibrate(s1, s2, klist)[source]

Generate the unimolecular reversible equilibrium reaction S1 <-> S2.

Parameters:
s1, s2Monomer or MonomerPattern

S1 and S2 in the above reaction.

klistlist of 2 Parameters or list of 2 numbers

Forward (S1 -> S2) and reverse rate constants (in that order). If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on the names and states of S1 and S2 and these parameters will be included at the end of the returned component list.

Returns:
componentsComponentSet

The generated components. Contains one reversible Rule and optionally two Parameters if klist was given as plain numbers.

Examples

Simple two-state equilibrium between A and B:

Model()
Monomer('A')
Monomer('B')
equilibrate(A(), B(), [1, 1])

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('A')
Monomer('A')
>>> Monomer('B')
Monomer('B')
>>> equilibrate(A(), B(), [1, 1]) 
ComponentSet([
 Rule('equilibrate_A_to_B', A() | B(), equilibrate_A_to_B_kf, equilibrate_A_to_B_kr),
 Parameter('equilibrate_A_to_B_kf', 1.0),
 Parameter('equilibrate_A_to_B_kr', 1.0),
 ])
pysb.macros.pore_bind(subunit, sp_site1, sp_site2, sc_site, size, cargo, c_site, klist)[source]

Generate rules to bind a monomer to a circular homomeric pore.

The pore structure is defined by the pore_species macro – subunit monomers bind to each other from sp_site1 to sp_site2 to form a closed ring. The binding reaction takes the form pore + cargo | pore:cargo.

Parameters:
subunitMonomer or MonomerPattern

Subunit of which the pore is composed.

sp_site1, sp_site2string

Names of the sites where one copy of subunit binds to the next.

sc_sitestring

Name of the site on subunit where it binds to the cargo cargo.

sizeinteger

Number of subunits in the pore at which binding will occur.

cargoMonomer or MonomerPattern

Cargo that binds to the pore complex.

c_sitestring

Name of the site on cargo where it binds to subunit.

klistlist of Parameters or numbers

List containing forward and reverse rate constants for the binding reaction (in that order). Rate constants should either be both Parameter objects or both numbers. If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on the subunit, the pore size and the cargo, and these parameters will be included at the end of the returned component list.

Examples

Specify that a cargo molecule can bind reversibly to a 3-membered pore:

Model()
Monomer('Unit', ['p1', 'p2', 'sc_site'])
Monomer('Cargo', ['c_site'])
pore_bind(Unit, 'p1', 'p2', 'sc_site', 3, 
          Cargo(), 'c_site', [1e-4, 1e-1, 1])

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('Unit', ['p1', 'p2', 'sc_site'])
Monomer('Unit', ['p1', 'p2', 'sc_site'])
>>> Monomer('Cargo', ['c_site'])
Monomer('Cargo', ['c_site'])
>>> pore_bind(Unit, 'p1', 'p2', 'sc_site', 3, 
...           Cargo(), 'c_site', [1e-4, 1e-1, 1]) 
ComponentSet([
 Rule('pore_bind_Unit_3_Cargo',
     MatchOnce(Unit(p1=3, p2=1, sc_site=None) %
         Unit(p1=1, p2=2, sc_site=None) %
         Unit(p1=2, p2=3, sc_site=None)) +
         Cargo(c_site=None) |
     MatchOnce(Unit(p1=3, p2=1, sc_site=4) %
         Unit(p1=1, p2=2, sc_site=None) %
         Unit(p1=2, p2=3, sc_site=None) %
         Cargo(c_site=4)),
     pore_bind_Unit_3_Cargo_kf, pore_bind_Unit_3_Cargo_kr),
 Parameter('pore_bind_Unit_3_Cargo_kf', 0.0001),
 Parameter('pore_bind_Unit_3_Cargo_kr', 0.1),
 ])
pysb.macros.pore_transport(subunit, sp_site1, sp_site2, sc_site, min_size, max_size, csource, c_site, cdest, ktable)[source]

Generate rules to transport cargo through a circular homomeric pore.

The pore structure is defined by the pore_species macro – subunit monomers bind to each other from sp_site1 to sp_site2 to form a closed ring. The transport reaction is modeled as a catalytic process of the form pore + csource | pore:csource >> pore + cdest

Parameters:
subunitMonomer or MonomerPattern

Subunit of which the pore is composed.

sp_site1, sp_site2string

Names of the sites where one copy of subunit binds to the next.

sc_sitestring

Name of the site on subunit where it binds to the cargo csource.

min_size, max_sizeinteger

Minimum and maximum number of subunits in the pore at which transport will occur.

csourceMonomer or MonomerPattern

Cargo “source”, i.e. the entity to be transported.

c_sitestring

Name of the site on csource where it binds to subunit.

cdestMonomer or MonomerPattern

Cargo “destination”, i.e. the resulting state after the transport event.

ktablelist of lists of Parameters or numbers

Table of forward, reverse and catalytic rate constants for the transport reactions. The outer list must be of length max_size - min_size + 1, and the inner lists must all be of length 3. In the outer list, the first element corresponds to the transport through the pore of size min_size and the last element to that of size max_size. Each inner list contains the forward, reverse and catalytic rate constants (in that order) for the corresponding transport reaction, and each of these pairs must comprise solely Parameter objects or solely numbers (never some of each). If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on the subunit, the pore size and the cargo, and these parameters will be included at the end of the returned component list.

Examples

Specify that a three-membered pore is capable of transporting cargo from the mitochondria to the cytoplasm:

Model()
Monomer('Unit', ['p1', 'p2', 'sc_site'])
Monomer('Cargo', ['c_site', 'loc'], {'loc':['mito', 'cyto']})
pore_transport(Unit, 'p1', 'p2', 'sc_site', 3, 3,
               Cargo(loc='mito'), 'c_site', Cargo(loc='cyto'),
               [[1e-4, 1e-1, 1]])

Generates two rules–one (reversible) binding rule and one transport rule–and the three associated parameters.

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('Unit', ['p1', 'p2', 'sc_site'])
Monomer('Unit', ['p1', 'p2', 'sc_site'])
>>> Monomer('Cargo', ['c_site', 'loc'], {'loc':['mito', 'cyto']})
Monomer('Cargo', ['c_site', 'loc'], {'loc': ['mito', 'cyto']})
>>> pore_transport(Unit, 'p1', 'p2', 'sc_site', 3, 3,
...                Cargo(loc='mito'), 'c_site', Cargo(loc='cyto'),
...                [[1e-4, 1e-1, 1]]) 
ComponentSet([
 Rule('pore_transport_complex_Unit_3_Cargomito',
     MatchOnce(Unit(p1=3, p2=1, sc_site=None) %
         Unit(p1=1, p2=2, sc_site=None) %
         Unit(p1=2, p2=3, sc_site=None)) +
         Cargo(c_site=None, loc='mito') |
     MatchOnce(Unit(p1=3, p2=1, sc_site=4) %
         Unit(p1=1, p2=2, sc_site=None) %
         Unit(p1=2, p2=3, sc_site=None) %
         Cargo(c_site=4, loc='mito')),
     pore_transport_complex_Unit_3_Cargomito_kf,
     pore_transport_complex_Unit_3_Cargomito_kr),
 Parameter('pore_transport_complex_Unit_3_Cargomito_kf', 0.0001),
 Parameter('pore_transport_complex_Unit_3_Cargomito_kr', 0.1),
 Rule('pore_transport_dissociate_Unit_3_Cargocyto',
     MatchOnce(Unit(p1=3, p2=1, sc_site=4) %
         Unit(p1=1, p2=2, sc_site=None) %
         Unit(p1=2, p2=3, sc_site=None) %
         Cargo(c_site=4, loc='mito')) >>
     MatchOnce(Unit(p1=3, p2=1, sc_site=None) %
         Unit(p1=1, p2=2, sc_site=None) %
         Unit(p1=2, p2=3, sc_site=None)) +
         Cargo(c_site=None, loc='cyto'),
     pore_transport_dissociate_Unit_3_Cargocyto_kc),
 Parameter('pore_transport_dissociate_Unit_3_Cargocyto_kc', 1.0),
 ])
pysb.macros.synthesize(species, ksynth)[source]

Generate a reaction which synthesizes a species.

Note that species must be “concrete”, i.e. the state of all sites in all of its monomers must be specified. No site may be left unmentioned.

Parameters:
speciesMonomer, MonomerPattern or ComplexPattern

The species to synthesize. If a Monomer, sites are considered as unbound and in their default state. If a pattern, must be concrete.

ksynthParameters or number

Synthesis rate. If a Parameter is passed, it will be used directly in the generated Rule. If a number is passed, a Parameter will be created with an automatically generated name based on the names and site states of the components of species and this parameter will be included at the end of the returned component list.

Returns:
componentsComponentSet

The generated components. Contains the unidirectional synthesis Rule and optionally a Parameter if ksynth was given as a number.

Examples

Synthesize A with site x unbound and site y in state ‘e’:

Model()
Monomer('A', ['x', 'y'], {'y': ['e', 'f']})
synthesize(A(x=None, y='e'), 1e-4)

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('A', ['x', 'y'], {'y': ['e', 'f']})
Monomer('A', ['x', 'y'], {'y': ['e', 'f']})
>>> synthesize(A(x=None, y='e'), 1e-4) 
ComponentSet([
 Rule('synthesize_Ae', None >> A(x=None, y='e'), synthesize_Ae_k),
 Parameter('synthesize_Ae_k', 0.0001),
 ])
pysb.macros.synthesize_degrade_table(table)[source]

Generate a table of synthesis and degradation reactions.

Given a list of species, calls the synthesize and degrade macros on each one. The species and the parameter values are passed as a list of lists (i.e. a table) with each inner list consisting of the species, forward and reverse rates (in that order).

Each species’ associated pair of rates may be either Parameters or numbers. If Parameters are passed, they will be used directly in the generated Rules. If numbers are passed, Parameters will be created with automatically generated names based on the names and states of the relevant species and these parameters will be included in the returned component list. To omit any individual reaction, pass None in place of the corresponding parameter.

Note that any species with a non-None synthesis rate must be “concrete”.

Parameters:
tablelist of lists

Table of species and rates, as described above.

Returns:
componentsComponentSet

The generated components. Contains the unidirectional synthesis and degradation Rules and optionally the Parameters for any rates given as numbers.

Examples

Specify synthesis and degradation reactions for A and B in a table:

Model()
Monomer('A', ['x', 'y'], {'y': ['e', 'f']})
Monomer('B', ['x'])
synthesize_degrade_table([[A(x=None, y='e'), 1e-4, 1e-6],
                          [B(),              None, 1e-7]])

Execution:

>>> Model() 
<Model '_interactive_' ...>
>>> Monomer('A', ['x', 'y'], {'y': ['e', 'f']})
Monomer('A', ['x', 'y'], {'y': ['e', 'f']})
>>> Monomer('B', ['x'])
Monomer('B', ['x'])
>>> synthesize_degrade_table([[A(x=None, y='e'), 1e-4, 1e-6],
...                           [B(),              None, 1e-7]]) 
ComponentSet([
    Rule('synthesize_Ae', None >> A(x=None, y='e'), synthesize_Ae_k),
    Parameter('synthesize_Ae_k', 0.0001),
    Rule('degrade_Ae', A(x=None, y='e') >> None, degrade_Ae_k),
    Parameter('degrade_Ae_k', 1e-06),
    Rule('degrade_B', B() >> None, degrade_B_k),
    Parameter('degrade_B_k', 1e-07),
    ])