Skip to content

Commit

Permalink
feat: add package extension for OrdinaryDiffEq (#150)
Browse files Browse the repository at this point in the history
  • Loading branch information
mehalter committed May 23, 2023
1 parent 856d763 commit 5d12a7a
Show file tree
Hide file tree
Showing 5 changed files with 34 additions and 5 deletions.
5 changes: 4 additions & 1 deletion Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -13,12 +13,14 @@ LabelledArrays = "2ee39098-c373-598a-b85f-a56591580800"
[weakdeps]
Catalyst = "479239e8-5488-4da2-87a7-35f2df7eef83"
ModelingToolkit = "961ee093-0014-501f-94e3-6117800e7a78"
OrdinaryDiffEq = "1dea7af3-3e70-54e6-95c3-0bf5283fa5ed"
Petri = "4259d249-1051-49fa-8328-3f8ab9391c33"

[extensions]
AlgebraicPetriPetriExt = "Petri"
AlgebraicPetriCatalystExt = "Catalyst"
AlgebraicPetriModelingToolkitExt = "ModelingToolkit"
AlgebraicPetriOrdinaryDiffEqExt = "OrdinaryDiffEq"
AlgebraicPetriPetriExt = "Petri"

[compat]
Catalyst = "13"
Expand All @@ -27,5 +29,6 @@ DataStructures = "0.18"
GeneralizedGenerated = "0.3"
LabelledArrays = "1"
ModelingToolkit = "8"
OrdinaryDiffEq = "6"
Petri = "1"
julia = "1.9"
9 changes: 9 additions & 0 deletions ext/AlgebraicPetriOrdinaryDiffEqExt.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,9 @@
module AlgebraicPetriOrdinaryDiffEqExt

using AlgebraicPetri
import OrdinaryDiffEq: ODEProblem

ODEProblem(m::AbstractPetriNet, u0, tspan, β) = ODEProblem(vectorfield(m), u0, tspan, β)
ODEProblem(m::AbstractPetriNet, tspan) = ODEProblem(vectorfield(m), concentrations(m), tspan, rates(m))

end
9 changes: 5 additions & 4 deletions src/AlgebraicPetri.jl
Original file line number Diff line number Diff line change
Expand Up @@ -13,7 +13,7 @@ export SchPetriNet, PetriNet, OpenPetriNetOb, AbstractPetriNet, ns, nt, ni, no,
Open, OpenPetriNet, OpenLabelledPetriNet, OpenReactionNet, OpenLabelledReactionNet,
OpenPetriNetOb, OpenLabelledPetriNetOb, OpenReactionNetOb, OpenLabelledReactionNetOb,
mca, flatten_labels,
AbstractPropertyPetriNet, sprop, tprop, sprops, tprops,
AbstractPropertyPetriNet, AbstractPropertyReactionNet, sprop, tprop, sprops, tprops,
SchPropertyPetriNet, SchPropertyLabelledPetriNet, SchPropertyReactionNet, SchPropertyLabelledReactionNet,
PropertyPetriNet, PropertyLabelledPetriNet, PropertyReactionNet, PropertyLabelledReactionNet,
OpenPropertyPetriNet, OpenPropertyLabelledPetriNet, OpenPropertyReactionNet, OpenPropertyLabelledReactionNet,
Expand Down Expand Up @@ -558,7 +558,7 @@ See Catlab.jl documentation for description of the @present syntax.
sname::Attr(S, Name)
end

@abstract_acset_type AbstractLabelledReactionNet <: AbstractPetriNet
@abstract_acset_type AbstractLabelledReactionNet <: AbstractReactionNet
@acset_type LabelledReactionNetUntyped(SchLabelledReactionNet, index=[:it, :is, :ot, :os]) <: AbstractLabelledReactionNet
const LabelledReactionNet{R,C} = LabelledReactionNetUntyped{R,C,Symbol}
const OpenLabelledReactionNetObUntyped, OpenLabelledReactionNetUntyped = OpenACSetTypes(LabelledReactionNetUntyped, :S)
Expand Down Expand Up @@ -736,7 +736,8 @@ See Catlab.jl documentation for description of the @present syntax.
sprop::Attr(S, Prop)
tprop::Attr(T, Prop)
end
@acset_type PropertyReactionNet(SchPropertyReactionNet, index=[:it, :is, :ot, :os]) <: AbstractPropertyPetriNet
@abstract_acset_type AbstractPropertyReactionNet <: AbstractPropertyPetriNet
@acset_type PropertyReactionNet(SchPropertyReactionNet, index=[:it, :is, :ot, :os]) <: AbstractPropertyReactionNet

const OpenPropertyReactionNetOb, OpenPropertyReactionNet = OpenACSetTypes(PropertyReactionNet, :S)
Open(p::PropertyReactionNet{R,C,T}, legs...) where {R,C,T} = OpenPropertyReactionNet{R,C,T}(p, map(l -> FinFunction(l, ns(p)), legs)...)
Expand All @@ -752,7 +753,7 @@ See Catlab.jl documentation for description of the @present syntax.
sprop::Attr(S, Prop)
tprop::Attr(T, Prop)
end
@acset_type PropertyLabelledReactionNetUntyped(SchPropertyLabelledReactionNet, index=[:it, :is, :ot, :os]) <: AbstractPropertyPetriNet
@acset_type PropertyLabelledReactionNetUntyped(SchPropertyLabelledReactionNet, index=[:it, :is, :ot, :os]) <: AbstractPropertyReactionNet
const PropertyLabelledReactionNet{R,C,T} = PropertyLabelledReactionNetUntyped{R,C,Symbol,T}

const OpenPropertyLabelledReactionNetObUntyped, OpenPropertyLabelledReactionNetUntyped = OpenACSetTypes(PropertyLabelledReactionNetUntyped, :S)
Expand Down
1 change: 1 addition & 0 deletions test/Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -5,6 +5,7 @@ Catalyst = "479239e8-5488-4da2-87a7-35f2df7eef83"
Catlab = "134e5e36-593f-5add-ad60-77f754baafbe"
LabelledArrays = "2ee39098-c373-598a-b85f-a56591580800"
ModelingToolkit = "961ee093-0014-501f-94e3-6117800e7a78"
OrdinaryDiffEq = "1dea7af3-3e70-54e6-95c3-0bf5283fa5ed"
Petri = "4259d249-1051-49fa-8328-3f8ab9391c33"
Tables = "bd369af6-aec1-5ad0-b16a-f7cc5008161c"
Test = "8dfed614-e22c-5e08-85e1-65c5234f0b40"
15 changes: 15 additions & 0 deletions test/ext/AlgebraicPetriOrdinaryDiffEqExt.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,15 @@
module TestAlgebraicPetriOrdinaryDiffEqExt

using Test
using AlgebraicPetri
using OrdinaryDiffEq

sir_petri = PetriNet(3, ((1, 2), (2, 2)), (2, 3))
sir_rxn = ReactionNet{Number,Int}([990, 10, 0], (0.001, ((1, 2) => (2, 2))), (0.25, (2 => 3)))

sol_petri = solve(ODEProblem(sir_petri, [990, 10, 0], (0,100.0), [0.001, 0.25]), Tsit5())
sol_rxn = solve(ODEProblem(sir_rxn, (0,100.0)), Tsit5())

@test sol_petri.u == sol_rxn.u

end

0 comments on commit 5d12a7a

Please sign in to comment.