AlgebraicRewriting.Rewrite.PBPO.PBPORuleType
  l    r

L ⟵ K ⟶ R tl ↓ ↓ tk <== tl, tk must be monic L' ⟵ K'

It is assumed we never want the typing/adherence match to be monic, but we can optionally restrict the match L → G to be monic.

We can attach application conditions to both the match morphism as well as the adherence morphism. Until morphism search under constraints becomes efficient, it's sometimes needed to just directly state the adherence morphism as a function of the match morphism.

AlgebraicRewriting.Rewrite.PBPO.canonMethod

Take a PBPO rule and put into normal form, i.e. where the lower square forms a pullback

See Prop 2.4 of "The PBPO graph transformation approach"

AlgebraicRewriting.Rewrite.PBPO.partial_abstractMethod

This construction addresses the following problem: ideally when we 'abstract' an ACSet from X to A->X, maps into X, say B->X, can be canonically pulled back to maps B->A which commute. However, A won't do here, because there may not even exist any maps B->A. If B has concrete attributes, then those cannot be sent to an AttrVar in A. Furthermore, if B has multiple 'references' to an AttrVar (two different edges, each with AttrVar(1), sent to two different edges with the same atttribute value in X), then there is no longer a canonical place to send AttrVar(1) to in A, as there is a distinct AttrVar for every single part+attr in X. So we need a construction which does two things to A->X, starting with a map B->X. 1.) replaces exactly the variables we need with concrete values in order to allow a map B->A, 2.) quotients variables in A so that there is exactly one choice for where to send attrvars in B such that the triangle commutes.

Starting with a map L -> G (where G has no AttrVars), we want the analogous map into a "partially abstracted" version of G that has concrete attributes replaced with AttrVars EXCEPT for those attributes which are mapped to by concrete attributes of L. Likewise, multiple occurences of the same variable in L correspond to AttrVars which should be merged in the partially-abstracted G.

For example, for a schema with a single Ob and Attr (where all combinatorial maps are just {1↦1, 2↦2}):

  • L = [AttrVar(1), :foo]

  • G = [:bar, :foo, :baz]

  • abs(G) = [AttrVar(1), AttrVar(2), AttrVar(3)]

  • expected result: [AttrVar(1), :foo, AttrVar(2)]

    L -> Partial_abs(G) ↓ ↑ G <- abs(G)

This function computes the top arrow of this diagram starting with the left arrow. The bottom arrow is computed by abstract_attributes and the right arrow by sub_vars. Furthermore, a map from Partial_abs(G) to G is provided.

This is the factorization system arising from a coreflective subcategory.

(see https://ncatlab.org/nlab/show/reflective+factorization+system and https://blog.algebraicjulia.org/post/2023/06/varacsets/)

AlgebraicRewriting.Rewrite.Utils.get_matchesMethod

PBPO matches consist of two morphisms. First, a match m: L → G, and secondly a typing G → L′. With attributes, it is not so simple because G has concrete values for attributes and L′ may have variables. Therefore, we actually change the typing to map out of A, an abstracted version of G (with its attributes replaced by variables). So we lift matches L->G to matches L->A, then search α∈Hom(A,L′).

In general, we want α to be uniquely determined by m, so by default α_unique is set to true.

 m

L⌟ ⟶ G || ↓ α L ⟶ L′ tl

 m

L ⟶ G tl ↓ ↘a ↑ (abs = partial abstraction. Note a is Labs in the code.) L′⟵ A α

The "strong match" condition we enforce is that: tl⁻¹(α(A)) = a⁻¹(A). This means we can deduce precisely what m is by looking at α.

AlgebraicRewriting.Rewrite.Utils.rewrite_match_mapsMethod
         r
      K ----> R
gₗ    u ↓ gᵣ ⌜ ↓ w

Gₗ <–– Gk ––> Gᵣ α ↓ ⌞ ↓ u' L′ <– K′ tₗ

For the adherence morphism α to be valid, it must satisfy a condition with m, tₗ. This is checked for matches provided by get_matches, so by default we do not check it.

L <–⌞• m ↓ ↓ G ⟵ Gk

See Lemma 7.2 of "TERMINATION OF GRAPH TRANSFORMATION SYSTEMS USING WEIGHTED SUBGRAPH COUNTING" by Overbeek and Endrullis (2023)

AlgebraicRewriting.Incremental.IncrementalHom.IncHomSetType

There are currently two types of IncHomSets. Common to both is a separation of the data required to specify which hom-set is being maintained (static) and the data structure required to be maintained as some target X is changing (runtime). We also distinguish between general static data for maintaining an incremental hom-set from further static data about various constraints on that process (constraints).

These are all currently designed to work with DenseParts ACSets but could be straightforwardly modified to work with MarkAsDeleted ACSets (which would be even more efficient.)

AlgebraicRewriting.Incremental.IncrementalHom.IncRuntimeType

Runtimes must have a state, key_vect, and key_dict.

state: The current state of the world, into which we are maintaining a hom-set.

key_vect::Vector{K}: A vector of keys into the ground source of truth, which stores the data of the morphisms. The structure of this source of truth and key type K depends on which kind of IncRuntime is being used. This may contain references to invalidated homs.

key_dict::Dict{K, Int}: Contains an entry for each current element of the hom-set. The values are used to index key_vect.

AlgebraicRewriting.Incremental.IncrementalHom.IncStaticType

A incremental hom-set has a fixed pattern, pattern (sometimes referred to by a variable, L). It is updated with respect to a class of changes to state (sometimes referred to by a variable X).

The state, X, can be changed via restriction to subobjects or by pushout with any of the predeclared monic addition morphisms, though we won't know the match morphism •→X until runtime.

                          additionⱼ
                        I >--------> R
                        ↓          ⌜ ↓
                        X ---------> X′

Cached overlaps data helps compute how the hom set updates with respect to a given addition being applied. The cache records in element i data for when we are updating the hom set after having applied addition i. It stores partial overlaps between L and Rᵢ which have some data that hits upon the new content added.

AlgebraicRewriting.Incremental.IncrementalSum.add_keys!Method

Add keys to a component of a IncSumHomSet. This implicitly adds composite keys for every combination of existing keys in the other components' hom sets.

Returns the composite keys which are added as an explicit list.

Base.lengthMethod

How many additions we've seen so far (including the initialization of the hom set). Could be confused with length(h.components) or length(keys(h))

AlgebraicRewriting.CategoricalAlgebra.CSets.check_pbMethod

Y i↘ f_ X → • g_ ↓ ⌟ ↓ f • → • g

Check whether (X, f,g) is the pullback of (f,g), up to isomorphism (i.e. the pullback of f and g produces (Y,π₁,π₂), where Y is isomorphic to X and i⋅f_ = π₁ & i⋅g_ = π₂.

AlgebraicRewriting.CategoricalAlgebra.CSets.sub_varsFunction

Given a value for each variable, create a morphism X → X′ which applies the substitution. We do this via pushout.

O –> X where C has AttrVars for merge equivalence classes ↓ and O has only AttrVars (sent to concrete values or eq classes C in the map to C.

subs and merge are dictionaries keyed by attrtype names

subs values are int-keyed dictionaries indicating binding, e.g. ; subs = (Weight = Dict(1 => 3.20, 5 => 2.32), ...)

merge values are vectors of vectors indicating equivalence classes, e.g. ; merge = (Weight = [[2,3], [4,6]], ...)

AlgebraicRewriting.CategoricalAlgebra.FinSets.pushout_complementMethod

Compute pushout complement of attributed C-sets, if possible.

The pushout complement is constructed pointwise from pushout complements of finite sets. If any of the pointwise identification conditions fail (in FinSet), this method will raise an error. If the dangling condition fails, the resulting C-set will be only partially defined. To check all these conditions in advance, use the function can_pushout_complement.

AlgebraicRewriting.Incremental.IncrementalHom.IncHomSetMethod

Automatically determine whether one creates an IncCC or an IncHom.

single keyword forces the pattern to be treated as a single CC, even if it's not one. Having any constraints will also force it to be treated as a single CC.

AlgebraicRewriting.Incremental.IncrementalHom.IncHomSetMethod

Initialize an Incremental hom set from a rule, using it's pattern L as the domain of the hom set. Any additions other than the map I->R can be passed in by additions, and the schema Presentation can be passed as S.

AlgebraicRewriting.Schedules.PolyModule

Mealy machines (augmented with monadic output) are a user-friendly format for specifying a behavior tree. Behavior trees in general are not finitely expressible, but we focus on trees which can be lazily generated by functions.

Although it is conceptually simple to think of a single set of "input doors" out "output doors" to enter/leave the Mealy machine, such that a Mealy machine has type A → B, we use Σᵢ Aᵢ → Σⱼ Bⱼ, where Aᵢ and Bⱼ are Julia types. This allows us to represent a Mealy machine with (Int + String)-many input doors, for example.

AlgebraicRewriting.Schedules.Poly.BTreeType

Lazily grown behavior tree induced by a Mealy machine. The future behavior is dictated by the inputs seen thus far (i.e. a vector of WireVals).

Each vertex is identified by a sequence of inputs and has a state of the Mealy machine associated with it.

Each nonempty sequence of inputs has a MealyRes associated with it.

AlgebraicRewriting.Schedules.Poly.MealyType

A function that maintains a state (initially s0) and has monadic output for some polynomial monad t.

The function f must be of type S × WireVal → S × (t ◁ WireVal)

                            (i.e. S × Wireval → MealyRes)
AlgebraicRewriting.Schedules.Poly.MealyResType

Output of a Mealy machine

newS - The new state of the Mealy machine mval - outputs along with their monadic values (e.g. probability weights) msg - A message reporting something about the computation

AlgebraicRewriting.Schedules.Poly.PMonadType

Quotes in this docstring and others taken from David Spivak: https://topos.site/blog/2023/09/powers-of-polynomial-monads/#exponentiating-monads

"A polynomial monad is a polynomial functor t with coherent maps η: y → t and μ: t ◁ t → t.

Polynomial monads can be thought of as offering compositional (possibly labeled) packages with some number of slots. The compositionality of this packaging says that (via the monad unit) we know how to package up a given element of any set, and that (via the monad multiplication) we can take a package of packages and simplify it to a single package."

We consider monads t of the form: t = Σ_{i ∈ t(1)} y^{t[i]}

I is the type of labels, e.g. probability weights.

AlgebraicRewriting.Schedules.Poly.SimulatorType

Equips a wiring diagram description of a simulator with mutable data structures (now behavior trees for each box, but possibly incremental homomorphism caches in the future). Requires that all the boxes of the WiringDiagram be convertable to BTrees.

AlgebraicRewriting.Schedules.Poly.apply_scheduleMethod

In theory applying a schedule should result in a list of ACSets associated with out ports and monad labels (e.g. probabilities), and if one were to want to recover the trajectory of the output one would have to use a Writer monad of some sort. For simplicity, the application of a schedule will simply return the trajectories themselves (monadic multiplication could in principle condense this output to the pure output).

AlgebraicRewriting.Schedules.Poly.joindistMethod

"The lotteries monad returns the set of packages labeled with a lottery (a natural number N and a probability distribution on it) and again containing N-many slots."

AlgebraicRewriting.Schedules.Queries.QueryType

Has an A input/output and a B input/output (by default, the B input can be changed to some other type if needed).

A  R ---------↖
↓  ↓          []

⌜–––-⌝ [] | Query | [agent subroutine] ⌞–––-⌟ [] ↓ ↓ ↓ [] A B ∅ [] ↘–––––-↗ Performs one action per element of Hom(B,X), optionally with some constraints. (i.e. sends you out along the B wire with agent Bₙ->X).

After you have done this for all Bₙ, then you exit the A port (you need to update the A->X map, and, if at any point the agent was deleted, then you exit a third door typed by 0).

A constraint optionally will be applied to (1) the A->W<-B cospan of old agent and purported new agent. (the new agent is the first argument to the constraint)

AlgebraicRewriting.Processes.find_depsMethod

For a concrete sequence of rewrites applications [a₁,a₂...aₙ], compute a poset on the set of applications which reflects their casual connections, where a < b mean that a must occur temporaly before b.

AlgebraicRewriting.CategoricalAlgebra.PartialMap.partial_map_classifier_universal_propertyMethod

A partial function is defined by the following span: m f A ↩ X → B

We compute ϕ(m,f): A ⟶ T(B) such that the following is a pullback square: f X ⟶ B m ↓ ↓ η(B) A ⟶ T(B) ϕ

Essentially, ϕ sends elements of A to the 'real' values in T(B) when A is in the subobject picked out by X. When A is 'deleted', it picks out the right element of the additional data added by T(B).

AlgebraicRewriting.CategoricalAlgebra.PartialMap.partial_map_functor_homMethod

Because the functorial embedding of objects keeps a copy of the original data, what to do with morphisms is just carry them along. Because our implementation adds all of the additional stuff afterwards, index-wise, we can use literally the same data for a morphism lifted from X⟶Y to T(X)⟶T(Y).

However, we still need to map the extra stuff in T(X) to the proper extra stuff in T(Y).

AlgebraicRewriting.CategoricalAlgebra.PartialMap.partial_map_functor_obMethod

A functor T, playing the role of Maybe in Set, but generalized to C-Sets.

When called on the terminal object, this produces the subobject classifier: See Mulry "Partial map classifiers and cartesian closed categories" (1994)

This function specifies what T does on objects. The key properties:

  1. for all X ∈ Ob(C), η(X):X⟶T(X) is monic. m f ϕ(m,f)
  2. for each span A ↩ X → B, there exists a unique morphism A ⟶ T(B) such that (m,f) is the pullback of ϕ(m,f),η(B))

Not only do we add an extra element to each component of the C-Set, but we need to consider the possibility that a component (with n outgoing morphisms) has any combination of the targets of those morphisms deleted (like the subobject classifier, there are different ways for something to be deleted).

For example, in Graph, an edge can be deleted that goes between any two vertices of the graph. We can't map all deleted edges to the same point in T(E) (if we're going to satisfy that desired property #2), so we need an extra edge in T(E) for every possibility (from V1 to V2, from V1 to V3, ..., from [Deleted] to V1, ..., from V2 to [Deleted], ... from [Deleted] to [Deleted]), where [Deleted] is our name for the extra element added to T(V).

                [src]     [tgt]

Thus, T(E) ≅ |E| + (|V|+1) × (|V|+1).

In general, T(X) ≅ |X| + ∏ₕ(|T(codom(h))|) for each outgoing morphism h::X⟶Y

  • the |X| corresponds to the 'real' elements of X
  • the second term corresponds to the possible ways an X can be deleted.
  • This recursive formula means we require the schema of the C-set to be acyclic otherwise the size is infinite (assumes schema is free).
AlgebraicRewriting.CategoricalAlgebra.StructuredCospans.composeH_Method
composeH_(r₁, r₂)

compose two rewrite rules horizontally (via pushouts) as shown below: L₁₋₍ₙ₋₁₎-> L <- Lₙ X₁ -> X <- X₂₋ₘ L₁₋₍ₙ₋₁₎ -> L +Lₙ X <- X₂₋ₘ ↑ λ ↑ ↑ ↑ ↑ χ ↑ ↑ ↑ ↑ I₁₋₍ₙ₋₁₎-> I <- Iₙ ∘h Y₁ -> Y <- Y₂₋ₘ = I₁₋₍ₙ₋₁₎ -> I +Iₙ Y <- Y₂₋ₘ ↓ ρ ↓ ↓ ↓ ↓ ζ ↓ ↓ ↓ ↓ R₁₋₍ₙ₋₁₎-> R <- Rₙ Z₁ -> Z <- Z₂₋ₘ R₁₋₍ₙ₋₁₎ -> R +Rₙ Z <- Z₂₋ₘ

AlgebraicRewriting.CategoricalAlgebra.StructuredCospans.composeV_Method
composeV_(r₁, r₂)

compose two rewrite rules vertically with pullbacks, as shown below: L₁₋ₙ -> L ↑ ↑ I₁₋ₙ -> I ↓ ↓ L₁₋ₙ -> L R₁₋ₙ -> R ↑ ↑ ∘v = I₁₋ₙ ×ᵣ₁₋ₙ Θ₁₋ₙ -> I ×ᵣ Θ Λ₁₋ₙ -> Λ ↓ ↓ ↑ ↑ Ω₁₋ₙ -> Ω Θ₁₋ₙ -> Θ ↓ ↓ Ω₁₋ₙ -> Ω

AlgebraicRewriting.CategoricalAlgebra.StructuredCospans.open_pushout_complementMethod

Initial data: 4 structured cospans + 3 cospan morphisms: μ, λ, ρ g G₁₋ₙ –> G ↑ l ↑ μ L₁₋ₙ –> L ↑ i ↑ λ I₁₋ₙ –> I ↓ r ↓ ρ R₁₋ₙ –> R

Computed data: 2 new structured cospans + 4 cospan morphisms: γ, η, ik, rh G₁₋ₙ G ↑ k ↑ γ ik I₁₋ₙ -> K₁₋ₙ –> K <– I ↓ h ↓ η rh R₁₋ₙ -> H₁₋ₙ –> H <– R In the context of the legs of a multicospan, the indices 1-n refer to the n legs of the cospan. In the context of a map of multicospans, there are 1-(n+1) maps, with the first one designating the map of the apexes. Hence it can make sense to have the elements: zip(legs, maps[2:end]) = [(legᵢ, mapᵢ), ...]

Catlab.CategoricalAlgebra.HomSearch.homomorphismsMethod

Find homomorphisms between structured cospans. These are constrained to be iso on the legs of the cospans. Solving this w/ homomorphism finding requires a dynamic acset, and the current hack will be replaced once those are available.

A homomorphism backend that uses SAT/SMT would also make this viable to do without hacking.

AlgebraicRewriting.Rewrite.Utils.rewrite_match_mapsMethod
rewrite_match_maps(r::Rule{:CoNeg}, m)

Apply a CoNegation rewrite rule (given as a span, L↩I->R) to a ACSet using a monic match morphism m which indicates where to apply the rewrite. l r L <- I -> R m ↓ ↓ ↓ G <- K -> H where K = ~L ∨ I

This works for any type that implements bi-Heyting logic operators ~ and ∨.

This is described here. Essentially, it is partway between DPO and SPO. Suppose the rule tries to delete two things, one of which satisfies the dangling condition, the other violates it. While DPO would fail to apply at all, and SPO would delete both things (cascading the deletion for the latter), co-negation rewriting would simply delete the item which can be deleted without cascading and ignore the other element.

It includes a quote which indicates that this method should work even when the match morphism isn't monic, if it satisfies the identification condition. Supporting this is not yet implemented.

Match morphisms which bind attribute variables are not monic, hence we this form of rewriting doesn't support VarACSets. Intuitively, it feels like this restriction could be relaxed.

AlgebraicRewriting.Schedules.Theories.ThTracedMonoidalWithBidiagonalsModule

No documentation found.

AlgebraicRewriting.Schedules.Theories.##4214 is of type Nothing.

Summary

struct Nothing

ThTracedMonoidalWithBidiagonals

      Ob::TYPE ⊣ []
      
      #= /juliateam/.julia/packages/GATlab/UL4sq/src/stdlib/theories/categories.jl:24 =#
      Hom(dom, codom)::TYPE ⊣ [dom::Ob, codom::Ob]
      
      #= /juliateam/.julia/packages/GATlab/UL4sq/src/stdlib/theories/categories.jl:33 =#
      compose(f, g)::Hom(a, c) ⊣ [a::Ob, b::Ob, c::Ob, f::Hom(a, b), g::Hom(b, c)]
      
      (assoc := compose(compose(f, g), h) == compose(f, compose(g, h))) ⊣ [a::Ob, b::Ob, c::Ob, d::Ob, f::Hom(a, b), g::Hom(b, c), h::Hom(c, d)]
      
      id(a)::Hom(a, a) ⊣ [a::Ob]
      
      (idl := compose(id(a), f) == f) ⊣ [a::Ob, b::Ob, f::Hom(a, b)]
      (idr := compose(f, id(b)) == f) ⊣ [a::Ob, b::Ob, f::Hom(a, b)]
      
      #= /juliateam/.julia/packages/Catlab/5M12F/src/theories/Monoidal.jl:28 =#
      otimes(A, B)::Ob ⊣ [A::Ob, B::Ob]
      otimes(f, g)::Hom(otimes(A, C), otimes(B, D)) ⊣ [A::Ob, B::Ob, C::Ob, D::Ob, f::Hom(A, B), g::Hom(C, D)]
      munit()::Ob ⊣ []
      (otimes(otimes(A, B), C) == otimes(A, otimes(B, C))) ⊣ [A::Ob, B::Ob, C::Ob]
      (otimes(munit(), A) == A) ⊣ [A::Ob]
      (otimes(A, munit()) == A) ⊣ [A::Ob]
      (otimes(otimes(f, g), h) == otimes(f, otimes(g, h))) ⊣ [A::Ob, B::Ob, C::Ob, X::Ob, Y::Ob, Z::Ob, f::Hom(A, X), g::Hom(B, Y), h::Hom(C, Z)]
      (otimes(id(munit()), f) == f) ⊣ [A::Ob, B::Ob, f::Hom(A, B)]
      (otimes(f, id(munit())) == f) ⊣ [A::Ob, B::Ob, f::Hom(A, B)]
      (compose(otimes(f, g), otimes(h, k)) == otimes(compose(f, h), compose(g, k))) ⊣ [A::Ob, B::Ob, C::Ob, X::Ob, Y::Ob, Z::Ob, f::Hom(A, B), h::Hom(B, C), g::Hom(X, Y), k::Hom(Y, Z)]
      (id(otimes(A, B)) == otimes(id(A), id(B))) ⊣ [A::Ob, B::Ob]
      
      braid(A, B)::Hom(otimes(A, B), otimes(B, A)) ⊣ [A::Ob, B::Ob]
      #= /juliateam/.julia/packages/Catlab/5M12F/src/theories/Monoidal.jl:87 =#
      (compose(braid(A, B), braid(B, A)) == id(otimes(A, B))) ⊣ [A::Ob, B::Ob]
      (braid(A, otimes(B, C)) == compose(otimes(braid(A, B), id(C)), otimes(id(B), braid(A, C)))) ⊣ [A::Ob, B::Ob, C::Ob]
      (braid(otimes(A, B), C) == compose(otimes(id(A), braid(B, C)), otimes(braid(A, C), id(B)))) ⊣ [A::Ob, B::Ob, C::Ob]
      (braid(A, munit()) == id(A)) ⊣ [A::Ob]
      (braid(munit(), A) == id(A)) ⊣ [A::Ob]
      (compose(otimes(f, g), braid(B, D)) == compose(braid(A, C), otimes(g, f))) ⊣ [A::Ob, B::Ob, C::Ob, D::Ob, f::Hom(A, B), g::Hom(C, D)]
      
      trace(X, A, B, f)::Hom(A, B) ⊣ [X::Ob, A::Ob, B::Ob, f::Hom(otimes(X, A), otimes(X, B))]
      
      mmerge(A)::Hom(otimes(A, A), A) ⊣ [A::Ob]
      #= /juliateam/.julia/packages/AlgebraicRewriting/RwqH4/src/schedules/Theories.jl:11 =#
      create(A)::Hom(munit(), A) ⊣ [A::Ob]
      #= /juliateam/.julia/packages/AlgebraicRewriting/RwqH4/src/schedules/Theories.jl:13 =#
AlgebraicRewriting.Rewrite.Constraints.QuantifierType

Quantified edge

e - which edge is filled in kind - Exists, Forall, or Exists! st - "such that", restrict the domain of quantification via a condition monic - restrict domain of quanitification to only monic matches

Catlab.Theories.:⊕Method

Combine two constraints disjunctively, sharing as much of the computation graph as possible.

Catlab.Theories.:⊗Method

Combine two constraints conjunctively, sharing as much of the computation graph as possible (i.e. pushout along the maximum common subgraph)

AlgebraicRewriting.CategoricalAlgebra.FinSets.id_conditionMethod

Check identification condition for pushout complement of finite sets.

The identification condition says that the functions do not map (1) both a deleted item and a preserved item in L to the same item in G or (2) two distinct deleted items to the same item. It is trivially satisfied for injective functions.

Returns pair of iterators of

(1) a nondeleted item that maps to a deleted item in G (2) a pair of distinct items in L that are deleted yet mapped to the same item in G.

AlgebraicRewriting.CategoricalAlgebra.FinSets.pushout_complementMethod

Compute a pushout complement of finite sets, if possible.

Given functions $l: I → L$ and $m: L → G$ to form a pushout square

l

L ← I m ↓ ↓k G ← K g

define the set $K := G / m(L / l(I))$ and take $g: K ↪ G$ to be the inclusion. Then the map $k: I → K$ is determined by the map $l⋅m: I → G$ from the requirement that the square commutes.

Pushout complements exist only if the identification condition is satisfied. An error will be raised if the pushout complement cannot be constructed. To check this in advance, use can_pushout_complement.

AlgebraicRewriting.Schedules.RuleApps.RuleAppType

Has the semantics of applying the rule to some match that is found (no guarantees on which one, which should be controlled by application conditions). If rewrite occurs, exit mode 1, else exit mode 2.

The agent is related to the L and R patterns of the rule. This can be done via a Span, or implicitly as a homomorphism into "I" of the rewrite rule, and alternatively just from the shape of the agent alone (if it is identical to I, take the id map, otherwise take the unique morphism into I).

AlgebraicRewriting.Rewrite.Utils.rewrite_match_mapsMethod
rewrite_match_maps(r::Rule{:DPO}, m)

Apply a DPO rewrite rule (given as a span, L<-I->R) to a ACSet using a match morphism m which indicates where to apply the rewrite. l r L <- I -> R m ↓ ↓ ↓ G <- K -> H

This works for any type that implements pushout_complement and pushout

AlgebraicRewriting.Incremental.Algorithms.compute_overlapsMethod

Find all partial maps from the pattern to the addition, with some restrictions:

  1. Something must be mapped into the newly added material.
  2. Anything in L incident to a part mapped onto newly added material must be mapped to newly added material
AlgebraicRewriting.Incremental.Algorithms.nac_overlapMethod

Goal: find overlaps between a negative application condition, N, and some current state of the world, X, such that we can efficiently check if something has been deleted generates a match by allowing it to satisfy the NAC. The data:

u: X ↢ X' (we are updating via some deletion) nac: L → N

X and X' are big. L and N are small. X ∖ X' is small.

Let η ↣ ~N be any subobject of ~N, the complement of N, which is our best approximation to N ∖ L.

We want matches η -> X that send everything in L to something in X'. But something in η must map to something not in X'. We look at all matches and then filter.

But we can do something more efficient than enumerating Hom(η, X), as we we need only generate matches into ~u (our best approximation to the part of X that got deleted, X ∖ X').

Thus, our process for finding overlaps requires only searching for morphisms between two things which are themselves pattern-sized.

AlgebraicRewriting.Incremental.IncrementalCCModule

An implementation for incremental hom sets where the pattern is a single connectect component (or otherwise has constraints that prevent the pattern from being broken up)

AlgebraicRewriting.Incremental.IncrementalCC.IncCCRuntimeType

This assumes the state of the world is changed in discrete updates.

match_vect: ground source of truth. It stores, for each update to state, what new matches were found at that time (+ still exist in the present state - hence the data must be a Dict{Int} rather than a vector, so that we can delete elements).

key_vect: way of indexing this list of dicts where the first element refers to an index of match_vect and the second element refers to the keys in the associated Dict. This gives us a single-integer value way of referring to every hom that has been seen, including the deleted ones. Thus the current # of matches is not length(key_vect).

key_dict: the relationship between match_vect and key_vect. There is an one element in the dictionary for every key (across all the dictionaries in match_vect). The keys of this dictionary are keys of match_vect, and the values are the keys of key_vect. So the current # of matches is length(key_dict), modulo any post-hoc constraints.

AlgebraicRewriting.Incremental.IncrementalHom.addition!Method

Add to matches based on an addition #i specified via a pushout (rmap, update)

             For all overlaps:  apex ↪ L
                                ↓      ⇣ (find all of these!)

Hom(L, Xold) => Hom(L, X) Rᵢ ⟶ X rmap ⌞ ↑ update Xold

However, we only want maps L → X where elements not in the image of L are all sent to X elements which outside of the image of rmap.

This is to avoid double-counting with a slightly bigger overlap which has already been calculated between L and Rᵢ.

Note, adding things can also invalidate old matches if there are negative application or dangling conditions.

Returns the 'keys' of the deleted matches and added matches.

AlgebraicRewriting.Incremental.IncrementalHom.deletion!Method

Delete / modify existing matches based on the target ACSet being permuted or reduced to a subobject. If a match touches upon something which is deleted, remove the match. Given X ↩ X′ we are updating Hom(L, X) => Hom(L, X′)

In the presence of negative application conditions / dangling condition, a deletion can also add new matches. When deletion is performed by DPO, we can know statically where to search for newly added morphisms. However, if we simply have an arbitrary deletion, this is harder.

A NAC has been deactivated if there exists a morphism N ⇾ X that sends all of L to the nondeleted part of X & some of N to the deleted portion. The deleted portion should be tiny compared to the not deleted portion, so it's pretty bad to search for all morphisms N->X and filter by those which have something which has something in the deleted portion. This means we ought to compute the overlaps on the fly in the non-DPO case.

The dpo flag signals that f was produced via pushout complement, such that any NAC can take advantage of its cached overlaps if it has them. If it is not nothing, it contains the morphism L ↢ I that was used in POC.

Returns the 'keys' of the deleted matches and added matches.

Base.lengthMethod

How many additions we've seen so far (including the initialization of the hom set). This could be confused with giving the total number of matches, which is instead length(keys(hset))

Base.push!Method

Consider a new addition (and compute its partial overlaps w/ the pattern)

AlgebraicRewriting.Schedules.Theories.ThTracedMonoidalWithBidiagonals.Meta.theoryConstant

No documentation found.

AlgebraicRewriting.Schedules.Theories.##4214 is of type Nothing.

Summary

struct Nothing

ThTracedMonoidalWithBidiagonals

      Ob::TYPE ⊣ []
      
      #= /juliateam/.julia/packages/GATlab/UL4sq/src/stdlib/theories/categories.jl:24 =#
      Hom(dom, codom)::TYPE ⊣ [dom::Ob, codom::Ob]
      
      #= /juliateam/.julia/packages/GATlab/UL4sq/src/stdlib/theories/categories.jl:33 =#
      compose(f, g)::Hom(a, c) ⊣ [a::Ob, b::Ob, c::Ob, f::Hom(a, b), g::Hom(b, c)]
      
      (assoc := compose(compose(f, g), h) == compose(f, compose(g, h))) ⊣ [a::Ob, b::Ob, c::Ob, d::Ob, f::Hom(a, b), g::Hom(b, c), h::Hom(c, d)]
      
      id(a)::Hom(a, a) ⊣ [a::Ob]
      
      (idl := compose(id(a), f) == f) ⊣ [a::Ob, b::Ob, f::Hom(a, b)]
      (idr := compose(f, id(b)) == f) ⊣ [a::Ob, b::Ob, f::Hom(a, b)]
      
      #= /juliateam/.julia/packages/Catlab/5M12F/src/theories/Monoidal.jl:28 =#
      otimes(A, B)::Ob ⊣ [A::Ob, B::Ob]
      otimes(f, g)::Hom(otimes(A, C), otimes(B, D)) ⊣ [A::Ob, B::Ob, C::Ob, D::Ob, f::Hom(A, B), g::Hom(C, D)]
      munit()::Ob ⊣ []
      (otimes(otimes(A, B), C) == otimes(A, otimes(B, C))) ⊣ [A::Ob, B::Ob, C::Ob]
      (otimes(munit(), A) == A) ⊣ [A::Ob]
      (otimes(A, munit()) == A) ⊣ [A::Ob]
      (otimes(otimes(f, g), h) == otimes(f, otimes(g, h))) ⊣ [A::Ob, B::Ob, C::Ob, X::Ob, Y::Ob, Z::Ob, f::Hom(A, X), g::Hom(B, Y), h::Hom(C, Z)]
      (otimes(id(munit()), f) == f) ⊣ [A::Ob, B::Ob, f::Hom(A, B)]
      (otimes(f, id(munit())) == f) ⊣ [A::Ob, B::Ob, f::Hom(A, B)]
      (compose(otimes(f, g), otimes(h, k)) == otimes(compose(f, h), compose(g, k))) ⊣ [A::Ob, B::Ob, C::Ob, X::Ob, Y::Ob, Z::Ob, f::Hom(A, B), h::Hom(B, C), g::Hom(X, Y), k::Hom(Y, Z)]
      (id(otimes(A, B)) == otimes(id(A), id(B))) ⊣ [A::Ob, B::Ob]
      
      braid(A, B)::Hom(otimes(A, B), otimes(B, A)) ⊣ [A::Ob, B::Ob]
      #= /juliateam/.julia/packages/Catlab/5M12F/src/theories/Monoidal.jl:87 =#
      (compose(braid(A, B), braid(B, A)) == id(otimes(A, B))) ⊣ [A::Ob, B::Ob]
      (braid(A, otimes(B, C)) == compose(otimes(braid(A, B), id(C)), otimes(id(B), braid(A, C)))) ⊣ [A::Ob, B::Ob, C::Ob]
      (braid(otimes(A, B), C) == compose(otimes(id(A), braid(B, C)), otimes(braid(A, C), id(B)))) ⊣ [A::Ob, B::Ob, C::Ob]
      (braid(A, munit()) == id(A)) ⊣ [A::Ob]
      (braid(munit(), A) == id(A)) ⊣ [A::Ob]
      (compose(otimes(f, g), braid(B, D)) == compose(braid(A, C), otimes(g, f))) ⊣ [A::Ob, B::Ob, C::Ob, D::Ob, f::Hom(A, B), g::Hom(C, D)]
      
      trace(X, A, B, f)::Hom(A, B) ⊣ [X::Ob, A::Ob, B::Ob, f::Hom(otimes(X, A), otimes(X, B))]
      
      mmerge(A)::Hom(otimes(A, A), A) ⊣ [A::Ob]
      #= /juliateam/.julia/packages/AlgebraicRewriting/RwqH4/src/schedules/Theories.jl:11 =#
      create(A)::Hom(munit(), A) ⊣ [A::Ob]
      #= /juliateam/.julia/packages/AlgebraicRewriting/RwqH4/src/schedules/Theories.jl:13 =#
AlgebraicRewriting.Rewrite.Utils.RuleType

Rewrite rules which are (usually) encoded as spans. The L structure encodes a pattern to be matched. The R morphism encodes a replacement pattern to be substituted in. They are related to each other by an interface I with maps: L ⟵ I ⟶ R

A semantics (DPO, SPO, CoNeg, or SqPO) must be chosen.

Control the match-finding process by specifying whether the match is intended to be monic or not, as well as an optional application condition(s)

AlgebraicRewriting.Rewrite.Utils.can_matchMethod

Returns nothing if the match is acceptable for rewriting according to the rule, otherwise returns the reason why it should be rejected

homsearch = if we know ahead of time that m was obtained m via automatic hom search, then we do not need to make certain checks

AlgebraicRewriting.Rewrite.Utils.get_matchesMethod

Get list of possible matches based on the constraints of the rule

This function has the same behavior as the generic get_matches, but it is more performant because we do not have to query all homomorphisms before finding a valid match, in case n=1.

AlgebraicRewriting.Schedules.Conditionals.ConditionalType

A primitive box in a NestedDWD which does not change the state but redirects it out of one of n wires.

It contains a function (A->X) -> ℝⁿ. This optionally depends on the internal state. This weights probability for n outports, conditional on the status of an ACSet. If the function just depends on X rather than the whole morphism, withagent is false. If the function does not depend on the internal state (assumed to be true iff initial state is nothing), then withstate is false.

The state and update function are by default trivial.

AlgebraicRewriting.Incremental.IncrementalConstraints.ACType

An application condition (w/r/t a pattern, L) is given by a morphism out of L. This can be given two semantics, based on the existence of a particular morphism making a triangle commute. We can furthermore impose monic constraints on that derived morphism.

AlgebraicRewriting.Incremental.IncrementalConstraints.IncConstraintsType

We may not be merely interested maintaining a set of matches Hom(L,X), but instead we care only about the monic morphisms, or morphisms that satisfy some positive/negative application conditions (PAC/NAC). In particular, NAC poses a new challenge: deletion can introduce new matches. There are various ways of dealing with this, one of which would require predeclaring deletion morphisms L ↢ I. However, that approach would only work for DPO. So a less efficient approach that uses only the data of X ↢ X′ and I→N is to search for all morphisms that send some part of N to a deleted part of X and all of L to the nondeleted part of X (this will find all of the new morphisms, only the new morphisms, but will possibly contain duplicates).

TODO add dangling field

AlgebraicRewriting.Incremental.IncrementalConstraints.NACType

A negative application condition L -> N means a match L -> X is invalid if there exists a commuting triangle.

We cache the subobjects of ~N (our best approximation to the things that are in N but not in L), as this can be taken advantage of in detecting new matches when we delete something.

dpo argument allows us to pass in a morphism I->L so that we can compute the staticly known overlaps beforehand, rather than at runtime. These are partial overlaps between N and L, which enumerate the possible ways a part of N can be sent to something that is deleted (part of L/I).

AlgebraicRewriting.Incremental.IncrementalConstraints.PACType

A positive application condition L -> P means a match L -> X is valid only if there does not exist a commuting triangle.

When we add something via some addition: O -> R, this could activate hitherto invalid matches which were blocked by a PAC. To detect these incrementally, we cache overlaps (indexed by possible additions) that store the ways in which the addition could intersect with P.

AlgebraicRewriting.Rewrite.Utils.rewrite_match_mapsMethod
rewrite_match_maps(r::Rule{:SqPO},m; pres::Union{Nothing, Presentation}=nothing)

Sesqui-pushout is just like DPO, except we use a final pullback complement instead of a pushout complement.

r.L  r.R

L <-⌞K -> R m ↓ ↓k ↓ r I <- • ->⌜O i o