Library Reference

AlgebraicPetri

AlgebraicPetri.SchReactionNetConstant

ACSet definition for a Petri net with rates on transitions and concentrations on states.

See Catlab.jl documentation for description of the @present syntax.

source
AlgebraicPetri.AbstractPetriNetType

Abstract type for C-sets that contain a petri net.

This type encompasses C-sets where the schema for graphs is a subcategory of C. This includes, for example, graphs, symmetric graphs, and reflexive graphs, but not half-edge graphs.

source
AlgebraicPetri.LabelledPetriNetMethod

LabelledPetriNet(n, ts::Vararg{Union{Pair,Tuple}})

Constructs a LabelledPetriNet object with state names as elements of n and labelled transitions described by ts. Transitions are given as transition_name=>((input_states)=>(output_states)).

A LabelledPetriNet modelling the SIR model with 3 states and 2 transitions can be constructed as follows:

LabelledPetriNet([:S, :I, :R], :inf=>((:S,:I)=>(:I,:I)), :rec=>(:I=>:R))
source
AlgebraicPetri.LabelledReactionNetMethod

LabelledReactionNet{R,C}(n, ts::Vararg{Union{Pair,Tuple}}) where {R,C}

Constructs a LabelledReactionNet object with labelled state concentrations as elements of n and labelled transitions described by ts. R is the data type used to store rates and C is the data type used to store concentrations.

Transitions are given as (t_name=>t_rate)=>((input_states)=>(output_states)).

A LabelledReactionNet modelling the SIR model with 3 states and 2 transitions, an initial population of 10 susceptible, 1 infected, 0 recovered and an infection rate of 0.5 and recovery rate of 0.1 can be constructed as follows:

LabelledReactionNet{Float64, Float64}([:S=>10,:I=>1,:R=>0], (:inf=>0.5)=>((:S,:I)=>(:I,:I)), (:rec=>0.1)=>(:I=>:R))
source
AlgebraicPetri.ReactionNetMethod

ReactionNet{R,C}(n, ts::Vararg{Union{Pair,Tuple}}) where {R,C}

Constructs a ReactionNet object with state concentrations as elements of n and transitions described by ts. R is the data type used to store rates and C is the data type used to store concentrations.

Transitions are given as transition_rate=>((input_states)=>(output_states)).

A ReactionNet modelling the SIR model with 3 states and 2 transitions, an initial population of 10 susceptible, 1 infected, 0 recovered and an infection rate of 0.5 and recovery rate of 0.1 can be constructed as follows:

ReactionNet{Float64, Float64}([10,1,0], 0.5=>((1,2)=>(2,2)), 0.1=>(2=>3))
source
AlgebraicPetri.TransitionMatricesType

TransitionMatrices

This data structure stores the transition matrix of an AbstractPetriNet object. This is primarily used for constructing the vectorfield representation of the Petri net.

source
Core.TypeMethod

(::AbstractPetriNet)(pn::AbstractPetriNet)

Cast one type of AbstractPetriNet to another. Any unrepresented parts will be nothing.

pn = PetriNet(3, (1,2)=>(2,2), 2=>3)
labelled_pn = LabelledPetriNet(pn)
source
Core.TypeMethod

(::AbstractPetriNet)(tm::TransitionMatrices)

Construct any AbstractPetriNet from a given transition matrice representation.

source
Core.TypeMethod

(::AbstractPropertyPetriNet)(pn::AbstractPetriNet, sprops, tprops)

Add properties to the states and transitions of a given Petri Net.

source
Core.TypeMethod

(::AbstractPetriNet)(n::Int, ts::Vararg{Union{Pair,Tuple}})

Constructs any AbstractPetriNet object with n states and transitions described by ts. Transitions are given as (input_states)=>(output_states).

A PetriNet modelling the SIR model with 3 states and 2 transitions can be constructed as follows:

PetriNet(3, (1,2)=>(2,2), 2=>3)
source
AlgebraicPetri.OpenMethod

Open(p::AbstractPetriNet, legs...)

Generates on OpenPetriNet with legs bundled as described by legs

source
AlgebraicPetri.OpenMethod

Open(p::AbstractPetriNet)

Converts a PetriNet to an OpenPetriNet where each state is exposed as a leg of the cospan. The OpenPetriNet can be composed over an undirected wiring diagram (see this blog post for a description of this compositional tooling)

source
AlgebraicPetri.add_input!Method

add_input!(p::AbstractPetriNet,t,s;kw...)

Add an input relationship to the Petri net between the transition t and species s.

Returns the ID of the input relationship

source
AlgebraicPetri.add_inputs!Method

add_inputs!(p::AbstractPetriNet,n,t,s;kw...)

Add input relationships to the Petri net between the transitions t and species s.

Returns the ID of the input relationship

source
AlgebraicPetri.add_output!Method

add_output!(p::AbstractPetriNet,t,s;kw...)

Add an output relationship to the Petri net between the transition t and species s.

Returns the ID of the input relationship

source
AlgebraicPetri.add_outputs!Method

add_outputs!(p::AbstractPetriNet,n,t,s;kw...)

Add output relationships to the Petri net between the transitions t and species s.

Returns the ID of the input relationship

source
AlgebraicPetri.add_species!Method

Add n species to the Petri net. Label and concentration can be provided depending on the kind of Petri net.

Returns the ID of the species

source
AlgebraicPetri.add_species!Method

Add a species to the Petri net. Label and concentration can be provided depending on the kind of Petri net.

Returns the ID of the species

source
AlgebraicPetri.add_transition!Method

Add a transition to the Petri net. Label and rate can be provided depending on the kind of Petri net.

Returns the ID of the transition

source
AlgebraicPetri.flatten_labelsMethod

flatten_labels(pn::AbstractPetriNet)

Takes a labelled Petri net or reaction net and flattens arbitrarily nested labels on the species and the transitions to a single symbol who's previously nested parts are separated by _.

source
AlgebraicPetri.snamesMethod

Names of species in a Petri net

Note that this returns indices if labels are not present in the PetriNet

source
AlgebraicPetri.tnamesMethod

Names of transitions in a Petri net

Note that this returns indices if labels are not present in the PetriNet

source
AlgebraicPetri.vectorfieldMethod

vectorfield(pn::AbstractPetriNet)

Generates a Julia function which calculates the vectorfield of the Petri net being simulated under the law of mass action.

The resulting function has a signature of the form f!(du, u, p, t) and can be passed to the DifferentialEquations.jl solver package.

source
AlgebraicPetri.vectorfield_exprMethod

vectorfield_expr(pn::AbstractPetriNet)

Generates a Julia expression which is then evaluated that calculates the vectorfield of the Petri net being simulated under the law of mass action.

The resulting function has a signature of the form f!(du, u, p, t) and can be passed to the DifferentialEquations.jl solver package.

source

AlgebraicPetri.BilayerNetworks

AlgebraicPetri.Epidemiology

AlgebraicPetri.Epidemiology.oapply_epiMethod

oapply_epi(ex, args...)

Generates a LabelledPetriNet under a composition pattern described by the undirected wiring diagram ex. This requires that the boxes in ex are only labelled with labels from the following set:

[:infection, :exposure, :illness, :recovery, :death]
source

AlgebraicPetri.ModelComparison

AlgebraicPetri.ModelComparison.compareMethod

compare Calculates all homomorphisms from a single Petri net (apex) to other Petri nets. This is returned as a Multispan with apex as the apex and the homomorphisms as the legs.

source
AlgebraicPetri.ModelComparison.petri_homomorphismsMethod

petri_homomorphisms Calculates all homomorphisms between two PetriNet models and returns an array of ACSetTransformations that describe these homomorphisms.

NOTE: This function does restrict to homomorphisms which preserve the transition signatures (number of input/output wires).

source

AlgebraicPetri.OpenTransitions

AlgebraicPetri.OpenTransitions.OpenTMethod

OpenT(p::AbstractPetriNet)

Converts a PetriNet to an OpenPetriNetT where each transition is exposed as a leg of the cospan. The OpenPetriNetT can be composed over an undirected wiring diagram.

source

AlgebraicPetri.TypedPetri

AlgebraicPetri.TypedPetri.add_reflexivesMethod

Modify a typed petri net to add "reflexive transitions". These are transitions which go from one species to itself, so they don't change the mass action semantics, but they are important for stratification.

The idea behind this is similar to the fact that the product of the graph –– with itself is

    / / /

/

One needs to add self-edges to each of the vertices in –– to get

–– | /| | / | | / | |/ | ––

Example:

add_reflexives(sird_typed, [[:strata], [:strata], [:strata], []], infectious_ontology)
source
AlgebraicPetri.TypedPetri.oapply_typedMethod

Takes in a labelled petri net and an undirected wiring diagram, where each of the boxes is labeled by a symbol that matches the label of a transition in the petri net. Then produces a petri net given by colimiting the transitions together, and returns the ACSetTransformation from that Petri net to the type system.

The tnames argument renames the transitions in the resulting typed Petri net.

source
AlgebraicPetri.TypedPetri.pairwise_id_typed_petriMethod

Make typed Petri net with 'identity' transformation between species pairs.

Assumes a single species type and a single transition type. For each pair of places, generate a transition consuming 1 input for each place and producing 1 output in each place.

source
AlgebraicPetri.TypedPetri.prim_petriMethod

Takes in a petri net and a transition in that petri net, constructs a petri net with just that transition, returns the acsettransformation from that into the type system.

source
AlgebraicPetri.TypedPetri.typed_productMethod

This takes the "typed product" of two typed Petri nets p1 and p2, which has

  • a species for every pair of a species in p1 and a species in p2 with the same type
  • a transition for every pair of a transitions in p1 and a species in p2 with the same type

This is the "workhorse" of stratification; this is what actually does the stratification, though you may have to "prepare" the petri nets first with add_reflexives

Returns a typed model, i.e. a map in Petri.

source

Package Extensions

Catalyst.jl

Catalyst.ReactionSystemType

ReactionSystem(pn::AbstractPetriNet)

Convert a general PetriNet to a ReactionSystem This conversion forgets any labels or rates provided, and converts all parameters and variables into symbols. It does preserve the ordering of transitions and states though (Transition 1 has a rate of k[1], state 1 has a concentration of S[1])

source

ModelngToolkit.jl

ModelingToolkit.ODESystemType

ODESystem(p::AbstractPetriNet; name=:PetriNet, kws...)

Convert a general PetriNet to an ODESystem This conversion forgets any labels or rates provided, and converts all parameters and variables into symbols. It does preserve the ordering of transitions and states though (Transition 1 has a rate of k[1], state 1 has a concentration of S[1])

source

ODESystem(bn::Union{AbstractLabelledBilayerNetwork,AbstractBilayerNetwork}; name=:BilayerNetwork, simplify = false)

Convert a general Bilayer Network to an ODESystem This conversion forgets any labels or rates provided, and converts all parameters and variables into symbols. It does preserve the ordering of transitions and states though (Transition 1 has a rate of k[1], state 1 has a concentration of S[1]). Note that symbolic simplification is not enable by default to preserve the bilayer structure. This is useful when one wants to convert symbolic expressions back to a bilayer network.

source

OrdinaryDiffEq.jl

SciMLBase.ODEProblemType

ODEProblem(m::AbstractPetriNet, u0, tspan, β)

Convert a general PetriNet to an ODEProblem This takes the initial conditions and rates as parameters and will ignore any initial conditions and rates embedded in the Petri Net if they exist.

source

ODEProblem(m::AbstractPetriNet, tspan)

Convert a general PetriNet to an ODEProblem This assumes the concentrations and weights are embedded in the Petri Net.

source