Computational graph for general feynman diagrams

API

FeynmanDiagram.ComputationalGraphs.FeynmanGraphType
mutable struct FeynmanGraph{F<:Number,W}

Computational graph representation of a (collection of) Feynman diagram(s). All Feynman diagrams should share the same set of external and internal vertices.

Members:

  • id::Int the unique hash id to identify the diagram
  • name::Symbol name of the diagram
  • orders::Vector{Int} orders associated with the Feynman graph, e.g., loop/derivative orders
  • properties::FeynmanProperties diagrammatic properties, e.g., the operator vertices and topology
  • subgraphs::Vector{FeynmanGraph{F,W}} vector of sub-diagrams
  • subgraph_factors::Vector{F} scalar multiplicative factors associated with each subdiagram
  • operator::DataType node operation (Sum, Prod, etc.)
  • weight::W weight of the diagram

Example:

julia> g1 = FeynmanGraph([]; vertices=[𝑓⁺(1),𝑓⁻(2)], external_indices=[1,2], external_legs=[true,true])
1:f⁺(1)|f⁻(2)=0.0

julia> g2 = FeynmanGraph([]; vertices=[𝑓⁺(3),𝑓⁻(4)], external_indices=[1,2], external_legs=[true,true])
2:f⁺(3)|f⁻(4)=0.0

julia> g = FeynmanGraph([g1,g2]; vertices=[𝑓⁺(1),𝑓⁻(2),𝑓⁺(3),𝑓⁻(4)], operator=ComputationalGraphs.Prod(), external_indices=[1,2,3,4], external_legs=[true,true,true,true])
3:f⁺(1)|f⁻(2)|f⁺(3)|f⁻(4)=0.0=Ⓧ (1,2)
source
FeynmanDiagram.ComputationalGraphs.GraphType
mutable struct Graph{F<:Number,W}

A representation of a computational graph, e.g., an expression tree, with type stable node data.

Members:

  • id::Int the unique hash id to identify the diagram
  • name::Symbol name of the diagram
  • orders::Vector{Int} orders associated with the graph, e.g., derivative orders
  • subgraphs::Vector{Graph{F,W}} vector of sub-diagrams
  • subgraph_factors::Vector{F} scalar multiplicative factors associated with each subgraph. Note that the subgraph factors may be manipulated algebraically. To associate a fixed multiplicative factor with this graph which carries some semantic meaning, use the factor argument instead.
  • operator::DataType node operation. Addition and multiplication are natively supported via operators Sum and Prod, respectively. Should be a concrete subtype of AbstractOperator.
  • weight::W the weight of this node
  • properties::Any extra information of Green's functions.

Example:

julia> g1 = Graph([])
1=0.0

julia> g2 = Graph([]; factor=2)
2⋅2.0=0.0

julia> g = Graph([g1, g2]; operator=ComputationalGraphs.Sum())
3=0.0=⨁ (1,2)
source
Base.:*Method
function Base.:*(c1, g2::Graph{F,W}) where {F,W}

Returns a graph representing the scalar multiplication c1*g2.

Arguments:

  • c1 scalar multiple
  • g2 Feynman graph
source
Base.:*Method
function Base.:*(c1, g2::Graph{F,W}) where {F,W}

Returns a graph representing the scalar multiplication c1*g2.

Arguments:

  • c1 scalar multiple
  • g2 computational graph
source
Base.:*Method
function Base.:*(g1::Graph{F,W}, c2) where {F,W}

Returns a graph representing the scalar multiplication g1*c2.

Arguments:

  • g1 Feynman graph
  • c2 scalar multiple
source
Base.:*Method
function Base.:*(g1::Graph{F,W}, c2) where {F,W}

Returns a graph representing the scalar multiplication g1*c2.

Arguments:

  • g1 computational graph
  • c2 scalar multiple
source
Base.:*Method
function Base.:*(g1::Graph{F,W}, g2::Graph{F,W}) where {F,W}

Returns a graph g1 * g2 representing the graph product between g1 and g2.

Arguments:

  • g1 first computational graph
  • g2 second computational graph
source
Base.:+Method
function Base.:+(g1::Graph{F,W}, g2::Graph{F,W}) where {F,W}

Returns a graph g1 + g2 representing the addition of two Feynman diagrams g2 with g1. Diagrams g1 and g2 must have the same diagram type, orders, and external vertices.

Arguments:

  • g1 first Feynman graph
  • g2 second Feynman graph
source
Base.:+Method
function Base.:+(g1::Graph{F,W}, g2::Graph{F,W}) where {F,W}

Returns a graph g1 + g2 representing the addition of g2 with g1. Graphs g1 and g2 must have the same orders.

Arguments:

  • g1 first computational graph
  • g2 second computational graph
source
Base.:-Method
function Base.:-(g1::Graph{F,W}, g2::Graph{F,W}) where {F,W}

Returns a graph g1 - g2 representing the subtraction of g2 from g1. Diagrams g1 and g2 must have the same diagram type, orders, and external vertices.

Arguments:

  • g1 first Feynman graph
  • g2 second Feynman graph
source
Base.:-Method
function Base.:-(g1::Graph{F,W}, g2::Graph{F,W}) where {F,W}

Returns a graph g1 - g2 representing the subtraction of g2 from g1. Graphs g1 and g2 must have the same orders.

Arguments:

  • g1 first computational graph
  • g2 second computational graph
source
Base.convertMethod
function Base.convert(::Type{G}, g::FeynmanGraph{F,W}) where {F,W,G<:Graph}

Converts a FeynmanGraph `g` into a Graph, discarding its Feynman properties.
After conversion, graph `g` is no longer guaranteed to be a valid (group of) Feynman diagram(s).

# Arguments:
- `g`  computational graph
source
Base.printMethod
print(io::IO, graph::AbstractGraph)

Write an un-decorated text representation of an AbstractGraph `graph` to the output stream `io`.
source
Base.showMethod
show(io::IO, graph::AbstractGraph; kwargs...)

Write a text representation of an AbstractGraph `graph` to the output stream `io`.
source
FeynmanDiagram.ComputationalGraphs.build_derivative_graphMethod
function build_derivative_graph(graphs::AbstractVector{G}, orders::NTuple{N,Int};
    nodes_id=nothing) where {G<:Graph,N}

Constructs a derivative graph using forward automatic differentiation with given graphs and derivative orders.

Arguments:

  • graphs: A vector of graphs.
  • orders::NTuple{N,Int}: A tuple indicating the orders of differentiation. N represents the number of independent variables to be differentiated.
  • nodes_id: Optional node IDs to indicate saving their derivative graph.

Returns:

  • A dictionary containing the dual derivative graphs for all indicated nodes.

If isnothing(nodes_id), indicated nodes include all leaf and root nodes. Otherwise, indicated nodes include all root nodes and other nodes from nodes_id.

source
FeynmanDiagram.ComputationalGraphs.burn_from_targetleaves!Method
function burn_from_targetleaves!(graphs::AbstractVector{G}, targetleaves_id::AbstractVector{Int}; verbose=0) where {G <: AbstractGraph}

Removes all nodes connected to the target leaves in-place via "Prod" operators.

Arguments:

  • graphs: A vector of graphs.
  • targetleaves_id::AbstractVector{Int}: Vector of target leafs' id.
  • verbose: Level of verbosity (default: 0).

Returns:

  • The id of a constant graph with a zero factor if any graph in graphs was completely burnt; otherwise, nothing.
source
FeynmanDiagram.ComputationalGraphs.count_expanded_operationMethod
function count_expanded_operation(g::G) where {G<:AbstractGraph}

Returns the total number of operations in the totally expanded version (without any parentheses in the mathematical expression) of the graph.

Arguments:

  • g::Graph: graph for which to find the total number of operations in its expanded version.
source
FeynmanDiagram.ComputationalGraphs.external_legsMethod
function external_legs(g::FeynmanGraph)

Returns a list of Boolean indices external_legs (::Vector{Bool}) indicating which external vertices of FeynmanGraph g have real legs (true: real leg, false: fake leg).

source
FeynmanDiagram.ComputationalGraphs.external_vertexMethod
function external_vertex(ops::OperatorProduct;
    name="", factor=one(_dtype.factor), weight=zero(_dtype.weight), operator=Sum())

Create a ExternalVertex-type FeynmanGraph from given OperatorProduct ops, including several quantum operators for an purely external vertex.

source
FeynmanDiagram.ComputationalGraphs.feynman_diagramMethod
function feynman_diagram(subgraphs::Vector{FeynmanGraph{F,W}}, topology::Vector{Vector{Int}}, perm_noleg::Union{Vector{Int},Nothing}=nothing;
    factor=one(_dtype.factor), weight=zero(_dtype.weight), name="", diagtype::DiagramType=GenericDiag()) where {F,W}

Create a FeynmanGraph representing feynman diagram from all subgraphs and topology (connections between vertices), where each ExternalVertex is given in vertices, while internal vertices are constructed with external legs of graphs in vertices, or simply OperatorProduct in vertices.

Arguments:

  • subgraphs::Vector{FeynmanGraph{F,W}} all subgraphs of the diagram. All external operators of subgraphs constitute all operators of the new diagram.
  • topology::Vector{Vector{Int}} topology of the diagram. Each Vector{Int} stores operators' index connected with each other (as a propagator).
  • perm_noleg::Union{Vector{Int},Nothing}=nothing permutation of all the nonleg external operators. By default, setting nothing means to use the default order from subgraphs.
  • factor::F overall scalar multiplicative factor for this diagram (e.g., permutation sign)
  • weight weight of the diagram
  • name name of the diagram
  • diagtype type of the diagram

Example:

julia> V = [𝑓⁺(1)𝑓⁻(2)𝜙(3), 𝑓⁺(4)𝑓⁻(5)𝜙(6), 𝑓⁺(7)𝑓⁻(8)𝜙(9)];
julia> g = feynman_diagram(interaction.(V), [[1, 5], [3, 9], [4, 8]], [3, 1, 2])
7:f⁺(1)f⁻(2)ϕ(3)|f⁺(4)f⁻(5)ϕ(6)|f⁺(7)f⁻(8)ϕ(9)=0.0=Ⓧ (1,2,3,4,5,6)

julia> g.subgraphs
6-element Vector{FeynmanGraph{Float64, Float64}}:
 1:f⁺(1)f⁻(2)ϕ(3)=0.0
 2:f⁺(4)f⁻(5)ϕ(6)=0.0
 3:f⁺(7)f⁻(8)ϕ(9)=0.0
 4:f⁺(1)|f⁻(5)⋅-1.0=0.0
 5:ϕ(3)|ϕ(9)=0.0
 6:f⁺(4)|f⁻(8)⋅-1.0=0.0
source
FeynmanDiagram.ComputationalGraphs.flatten_all_chains!Method
function flatten_all_chains!(g::AbstractGraph; verbose=0)

F Flattens all nodes representing trivial unary chains in-place in the given graph g.

Arguments:

  • graphs: The graph to be processed.
  • verbose: Level of verbosity (default: 0).

Returns:

  • The mutated graph g with all chains flattened.
source
FeynmanDiagram.ComputationalGraphs.flatten_all_chains!Method
function flatten_all_chains!(graphs::Union{Tuple,AbstractVector{<:AbstractGraph}}; verbose=0)

Flattens all nodes representing trivial unary chains in-place in the given graphs.

Arguments:

  • graphs: A collection of graphs to be processed.
  • verbose: Level of verbosity (default: 0).

Returns:

  • The mutated collection graphs with all chains in each graph flattened.
source
FeynmanDiagram.ComputationalGraphs.flatten_chains!Method
function flatten_chains!(g::AbstractGraph)

Recursively flattens chains of subgraphs within the given graph `g` by merging certain trivial unary subgraphs 
into their parent graphs in the in-place form.

Acts only on subgraphs of `g` with the following structure: 𝓞 --- 𝓞' --- ⋯ --- 𝓞'' ⋯ (!),
where the stop-case (!) represents a leaf, a non-trivial unary operator 𝓞'''(g) != g, or a non-unary operation.

Arguments:

  • g::AbstractGraph: graph to be modified
source
FeynmanDiagram.ComputationalGraphs.flatten_chainsMethod
function flatten_chains(g::AbstractGraph) 

Recursively flattens chains of subgraphs within a given graph `g` by merging certain trivial unary subgraphs into their parent graphs,
This function returns a new graph with flatten chains, derived from the input graph `g` remaining unchanged.

Arguments:

  • g::AbstractGraph: graph to be modified
source
FeynmanDiagram.ComputationalGraphs.flatten_prod!Method
flatten_prod!(graph::G; map::Dict{Int,G}=Dict{Int,G}()) where {G<:AbstractGraph}

Recursively merge multi-product sub-branches within the given graph `g  by merging  product subgraphs 
into their parent product graphs in the in-place form.

Arguments:

  • g::AbstractGraph: graph to be modified
  • map::Dict{Int,G}=Dict{Int,G}(): A dictionary that maps the id of an original node with its corresponding new node after transformation.

In recursive transform, nodes can be visited several times by different parents. This map keeps track of those visited, and reuse those transformed sub-branches instead of recreating them.

source
FeynmanDiagram.ComputationalGraphs.flatten_sum!Method
flatten_sum!(graph::G; map::Dict{Int,G}=Dict{Int,G}()) where {G<:AbstractGraph}

Recursively merge multi-product sub-branches within the given graph `g  by merging  sum subgraphs 
into their parent sum graphs in the in-place form.

Arguments:

  • g::AbstractGraph: graph to be modified
  • map::Dict{Int,G}=Dict{Int,G}(): A dictionary that maps the id of an original node with its corresponding new node after transformation.

In recursive transform, nodes can be visited several times by different parents. This map keeps track of those visited, and reuse those transformed sub-branches instead of recreating them.

source
FeynmanDiagram.ComputationalGraphs.forwardADMethod
function forwardAD(diag::Graph{F,W}, ID::Int) where {F,W}

Given a  graph G and an id of graph g, calculate the derivative d G / d g by forward propagation algorithm.

Arguments:

  • diag::Graph{F,W}: Graph G to be differentiated
  • ID::Int: ID of the graph g
source
FeynmanDiagram.ComputationalGraphs.forwardAD_root!Method
function forwardAD_root!(graphs::AbstractVector{G}, idx::Int=1,
    dual::Dict{Tuple{Int,NTuple{N,Bool}},G}=Dict{Tuple{Int,Tuple{Bool}},G}()) where {G<:Graph,N}

Computes the forward automatic differentiation (AD) of the given graphs beginning from the roots.

Arguments:

  • graphs: A vector of graphs.
  • idx: Index for differentiation (default: 1).
  • dual: A dictionary that holds the result of differentiation.

Returns:

  • The dual dictionary populated with all differentiated graphs, including the intermediate AD.
source
FeynmanDiagram.ComputationalGraphs.groupMethod
function group(gv::AbstractVector{G}, indices::Vector{Int}) where {G<:FeynmanGraph}

Group the graphs in gv by the external operators at the indices indices. Return a dictionary of Vector{OperatorProduct} to GraphVector.

Example

julia> p1 = propagator(𝑓⁺(1)𝑓⁻(2));

julia> p2 = propagator(𝑓⁺(1)𝑓⁻(3));

julia> p3 = propagator(𝑓⁺(2)𝑓⁻(3));

julia> gv = [p1, p2, p3];

julia> ComputationalGraphs.group(gv, [1, 2])
Dict{Vector{OperatorProduct}, Vector{FeynmanGraph{Float64, Float64}}} with 3 entries:
  [f⁻(2), f⁺(1)] => [1:f⁺(1)|f⁻(2)⋅-1.0=0.0]
  [f⁻(3), f⁺(1)] => [2:f⁺(1)|f⁻(3)⋅-1.0=0.0]
  [f⁻(3), f⁺(2)] => [3:f⁺(2)|f⁻(3)⋅-1.0=0.0]

julia> ComputationalGraphs.group(gv, [1, ])
Dict{Vector{OperatorProduct}, Vector{FeynmanGraph{Float64, Float64}}} with 2 entries:
  [f⁻(3)] => [2:f⁺(1)|f⁻(3)⋅-1.0=0.0, 3:f⁺(2)|f⁻(3)⋅-1.0=0.0]
  [f⁻(2)] => [1:f⁺(1)|f⁻(2)⋅-1.0=0.0]

julia> ComputationalGraphs.group(gv, [2, ])
Dict{Vector{OperatorProduct}, Vector{FeynmanGraph{Float64, Float64}}} with 2 entries:
  [f⁺(2)] => [3:f⁺(2)|f⁻(3)⋅-1.0=0.0]
  [f⁺(1)] => [1:f⁺(1)|f⁻(2)⋅-1.0=0.0, 2:f⁺(1)|f⁻(3)⋅-1.0=0.0]
source
FeynmanDiagram.ComputationalGraphs.has_zero_subfactorsMethod
function has_zero_subfactors(g::AbstractGraph, operator_type::Type{<:AbstractOperator})

Determines whether the graph `g` has only zero-valued subgraph factors based on the specified operator type.
This function does not recurse through the subgraphs of `g`, so it only checks the immediate subgraph factors.
If `g` is a leaf (i.e., has no subgraphs), the function returns `false` by convention.

The behavior of the function depends on the operator type:
- `Sum`: Checks if all subgraph factors are zero.
- `Prod`: Checks if any subgraph factor is zero.
- `Power{N}`: Checks if the first subgraph factor is zero.
- Other `AbstractOperator`: Defaults to return `false`.

Arguments:

  • g::AbstractGraph: graph to be analyzed
  • operator: the operator used in graph g
source
FeynmanDiagram.ComputationalGraphs.interactionMethod
function interaction(ops::OperatorProduct; name="", reorder::Union{Function,Nothing}=nothing,
    factor=one(_dtype.factor), weight=zero(_dtype.weight), operator=Sum())

Create a Interaction-type FeynmanGraph from given OperatorProduct ops, including several quantum operators for a vertex. One can call a reorder function for the operators ordering.

source
FeynmanDiagram.ComputationalGraphs.linear_combinationMethod
function linear_combination(graphs::Vector{FeynmanGraph{F,W}}, constants::AbstractVector=ones(F, length(graphs))) where {F,W}

Given a vector 𝐠 of graphs each with the same type and external/internal vertices and an equally-sized vector 𝐜 of constants, returns a new graph representing the linear combination (𝐜 ⋅ 𝐠). The function identifies unique graphs from the input graphs and sums their associated constants. All input Graphs must have the same diagram type, orders, and external vertices.

Arguments:

  • graphs vector of input FeymanGraphs
  • constants vector of scalar multiples (defaults to ones(F, length(graphs))).

Returns:

  • A new FeynmanGraph{F,W} object representing the linear combination of the unique input graphs weighted by the constants,

where duplicate graphs in the input graphs are combined by summing their associated constants.

Example:

Given graphs `g1`, `g2`, `g1` and constants `c1`, `c2`, `c3`, the function computes `(c1+c3)*g1 + c2*g2`.
source
FeynmanDiagram.ComputationalGraphs.linear_combinationMethod
function linear_combination(graphs::Vector{Graph{F,W}}, constants::AbstractVector=ones(F, length(graphs))) where {F,W}

Given a vector 𝐠 of graphs and an equally-sized vector 𝐜 of constants, returns a new graph representing the linear combination (𝐜 ⋅ 𝐠). The function identifies unique graphs from the input graphs and sums their associated constants. All input graphs must have the same orders.

Arguments:

  • graphs vector of computational graphs
  • constants vector of scalar multiples (defaults to ones(F, length(graphs))).

Returns:

  • A new Graph{F,W} object representing the linear combination of the unique input graphs weighted by the constants,

where duplicate graphs in the input graphs are combined by summing their associated constants.

Example:

Given graphs `g1`, `g2`, `g1` and constants `c1`, `c2`, `c3`, the function computes `(c1+c3)*g1 + c2*g2`.
source
FeynmanDiagram.ComputationalGraphs.linear_combinationMethod
function linear_combination(g1::FeynmanGraph{F,W}, g2::FeynmanGraph{F,W}, c1, c2) where {F,W}

Returns a graph representing the linear combination c1*g1 + c2*g2. If g1 == g2, it will return a graph representing (c1+c2)*g1 Feynman Graphs g1 and g2 must have the same diagram type, orders, and external vertices.

Arguments:

  • g1 first Feynman graph
  • g2 second Feynman graph
  • c1: first scalar multiple (defaults to 1).
  • c2: second scalar multiple (defaults to 1).
source
FeynmanDiagram.ComputationalGraphs.linear_combinationMethod
function linear_combination(g1::Graph{F,W}, g2::Graph{F,W}, c1, c2) where {F,W}

Returns a graph representing the linear combination c1*g1 + c2*g2. If g1 == g2, it will return a graph representing (c1+c2)*g1. Graphs g1 and g2 must have the same orders.

Arguments:

  • g1 first computational graph
  • g2 second computational graph
  • c1 first scalar multiple
  • c2 second scalar multiple
source
FeynmanDiagram.ComputationalGraphs.mask_zero_subgraph_factorsMethod
function mask_zero_subgraph_factors(operator::Type{<:AbstractOperator}, subg_fac::Vector{F}) where {F}

Returns a list of indices that should be considered when performing the operation (e.g., Sum, Prod, Power), effectively masking out zero values as appropriate.

The behavior of the function depends on the operator type:
- `Sum`: Returns all indices that are not equal to zero.
- `Prod`: Returns the index of the first zero value, or all indices if none are found.
- `Power`: Returns `[1]`, or error if the power is negative.
- Other `AbstractOperator`: Defaults to return all indices.
source
FeynmanDiagram.ComputationalGraphs.merge_all_linear_combinations!Method
function merge_all_linear_combinations!(g::AbstractGraph; verbose=0)

Merges all nodes representing a linear combination of a non-unique list of subgraphs in-place in the given graph `g`.

Arguments:

  • g: An AbstractGraph.
  • verbose: Level of verbosity (default: 0).

Returns:

  • Optimized graph.

source
FeynmanDiagram.ComputationalGraphs.merge_all_linear_combinations!Method
function merge_all_linear_combinations!(graphs::Union{Tuple,AbstractVector{<:AbstractGraph}}; verbose=0)

Merges all nodes representing a linear combination of a non-unique list of subgraphs in-place in the given graphs.

Arguments:

  • graphs: A collection of graphs to be processed.
  • verbose: Level of verbosity (default: 0).

Returns:

  • Optimized graphs.

source
FeynmanDiagram.ComputationalGraphs.merge_all_multi_products!Method
function merge_all_multi_products!(g::Graph; verbose=0)

Merges all nodes representing a multi product of a non-unique list of subgraphs in-place in the given graph `g`.

Arguments:

  • g::Graph: A Graph.
  • verbose: Level of verbosity (default: 0).

Returns:

  • Optimized graph.

source
FeynmanDiagram.ComputationalGraphs.merge_all_multi_products!Method
function merge_all_multi_products!(graphs::Union{Tuple,AbstractVector{<:Graph}}; verbose=0)

Merges all nodes representing a multi product of a non-unique list of subgraphs in-place in the given graphs.

Arguments:

  • graphs: A collection of graphs to be processed.
  • verbose: Level of verbosity (default: 0).

Returns:

  • Optimized graphs.

source
FeynmanDiagram.ComputationalGraphs.merge_linear_combination!Method
function merge_linear_combination!(g::AbstractGraph)

Modifies a computational graph `g` by factorizing multiplicative prefactors, e.g.,
3*g1 + 5*g2 + 7*g1 + 9*g2 ↦ 10*g1 + 14*g2 = linear_combination(g1, g2, 10, 14).
Returns a linear combination of unique subgraphs and their total prefactors. 
Does nothing if the graph `g` does not represent a Sum operation.

Arguments:

  • g::AbstractGraph: graph to be modified
source
FeynmanDiagram.ComputationalGraphs.merge_linear_combinationMethod
function merge_linear_combination(g::AbstractGraph)

Returns a copy of computational graph `g` with multiplicative prefactors factorized,
e.g., 3*g1 + 5*g2 + 7*g1 + 9*g2 ↦ 10*g1 + 14*g2 = linear_combination(g1, g2, 10, 14).
Returns a linear combination of unique subgraphs and their total prefactors. 
Does nothing if the graph `g` does not represent a Sum operation.

Arguments:

  • g::AbstractGraph: graph to be modified
source
FeynmanDiagram.ComputationalGraphs.merge_multi_product!Method
function merge_multi_product!(g::Graph{F,W}) where {F,W}

Merge multiple products within a computational graph `g` if they share the same operator (`Prod`).
If `g.operator == Prod`, this function will merge `N` identical subgraphs into a single subgraph with a power operator `Power(N)`. 
The function ensures each unique subgraph is counted and merged appropriately, preserving any distinct subgraph_factors associated with them.

Arguments:

  • g::Graph: graph to be modified

Returns

  • A merged computational graph with potentially fewer subgraphs if there were repeating subgraphs with the Prod operator. If the input graph's operator isn't Prod, the function returns the input graph unchanged.
source
FeynmanDiagram.ComputationalGraphs.merge_multi_productMethod
function merge_multi_product(g::Graph{F,W}) where {F,W}

Returns a copy of computational graph `g` with multiple products merged if they share the same operator (`Prod`).
If `g.operator == Prod`, this function will merge `N` identical subgraphs into a single subgraph with a power operator `Power(N)`. 
The function ensures each unique subgraph is counted and merged appropriately, preserving any distinct subgraph_factors associated with them.

Arguments:

  • g::Graph: graph to be modified

Returns

  • A merged computational graph with potentially fewer subgraphs if there were repeating subgraphs with the Prod operator. If the input graph's operator isn't Prod, the function returns the input graph unchanged.
source
FeynmanDiagram.ComputationalGraphs.multi_productMethod
multi_product(graphs::Vector{Graph{F,W}}, constants::AbstractVector=ones(F, length(graphs))) where {F,W,C}

Construct a product graph from multiple input graphs, where each graph can be weighted by a constant. For graphs that are repeated more than once, it adds a power operator to the subgraph to represent the repetition. Moreover, it optimizes any trivial unary operators in the resulting product graph.

Arguments:

  • graphs::Vector{Graph{F,W}}: A vector of input graphs to be multiplied.
  • constants::AbstractVector: A vector of scalar multiples. If not provided, it defaults to a vector of ones of the same length as graphs.

Returns:

  • A new product graph with the unique subgraphs (or powered versions thereof) and the associated constants as subgraph factors.

Example:

Given graphs `g1`, `g2`, `g1` and constants `c1`, `c2`, `c3`, the function computes `(c1*c3)*(g1)^2 * c2*g2`.
source
FeynmanDiagram.ComputationalGraphs.multi_productMethod
function multi_product(g1::Graph{F,W}, g2::Graph{F,W}, c1=F(1), c2=F(1)) where {F,W,C}

Returns a graph representing the multi product c1*g1 * c2*g2. If g1 == g2, it will return a graph representing c1*c2 * (g1)^2 with Power(2) operator.

Arguments:

  • g1: first computational graph
  • g2: second computational graph
  • c1: first scalar multiple (defaults to 1).
  • c2: second scalar multiple (defaults to 1).
source
FeynmanDiagram.ComputationalGraphs.open_parenthesis!Method
open_parenthesis!(graph::G; map::Dict{Int,G}=Dict{Int,G}()) where {G<:AbstractGraph}

Recursively open parenthesis of subgraphs within the given graph `g`with in place form.  The graph eventually becomes 
a single Sum root node with multiple subgraphs that represents multi-product of nodes (not flattened).

Arguments:

  • g::AbstractGraph: graph to be modified
  • map::Dict{Int,G}=Dict{Int,G}(): A dictionary that maps the id of an original node with its corresponding new node after transformation.

In recursive transform, nodes can be visited several times by different parents. This map keeps track of those visited, and reuse those transformed sub-branches instead of recreating them. parents

source
FeynmanDiagram.ComputationalGraphs.optimize!Method
function optimize!(graphs::Union{Tuple,AbstractVector{<:AbstractGraph}}; level=0, verbose=0, normalize=nothing)

In-place optimization of given `graphs`. Removes duplicated leaves, flattens chains, 
merges linear combinations, and removes zero-valued subgraphs. When `level > 0`, also removes duplicated intermediate nodes.

Arguments:

  • graphs: A tuple or vector of graphs.
  • level: Optimization level (default: 0). A value greater than 0 triggers more extensive but slower optimization processes, such as removing duplicated intermediate nodes.
  • verbose: Level of verbosity (default: 0).
  • normalize: Optional function to normalize the graphs (default: nothing).

Returns

  • Returns the optimized graphs. If the input graphs is empty, it returns nothing.
source
FeynmanDiagram.ComputationalGraphs.optimizeMethod
function optimize(graphs::Union{Tuple,AbstractVector{<:AbstractGraph}}; level=0, verbose=0, normalize=nothing)

Optimizes a copy of given `graphs`. Removes duplicated nodes (when `level > 0`) or leaves, flattens chains, 
merges linear combinations, and removing zero-valued subgraphs.

Arguments:

  • graphs: A tuple or vector of graphs.
  • level: Optimization level (default: 0). A value greater than 0 triggers more extensive but slower optimization processes, such as removing duplicated nodes.
  • verbose: Level of verbosity (default: 0).
  • normalize: Optional function to normalize the graphs (default: nothing).

Returns:

  • A tuple/vector of optimized graphs.
source
FeynmanDiagram.ComputationalGraphs.plot_treeMethod
function plot_tree(graph::AbstractGraph; verbose = 0, maxdepth = 6)

Visualize the computational graph as a tree using ete3 python package

#Arguments

  • graph::AbstractGraph : the computational graph struct to visualize
  • verbose=0 : the amount of information to show
  • maxdepth=6 : deepest level of the computational graph to show
source
FeynmanDiagram.ComputationalGraphs.propagatorMethod
function propagator(ops::Union{OperatorProduct,Vector{QuantumOperator}};
    name="", factor=one(_dtype.factor), weight=zero(_dtype.weight), operator=Sum())

Create a Propagator-type FeynmanGraph from given OperatorProduct or Vector{QuantumOperator} ops, including two quantum operators.

source
FeynmanDiagram.ComputationalGraphs.relabel!Method
function relabel!(g::FeynmanGraph, map::Dict{Int,Int})

Relabels the quantum operators in `g` and its subgraphs according to `map`.
For example, `map = {1=>2, 3=>2}`` will find all quantum operators with labels 1 and 3, and then map them to 2.

Arguments:

  • g::FeynmanGraph: graph to be modified
  • map: mapping from old labels to the new ones
source
FeynmanDiagram.ComputationalGraphs.relabelMethod
function relabel(g::FeynmanGraph, map::Dict{Int,Int})

Returns a copy of `g` with quantum operators in `g` and its subgraphs relabeled according to `map`.
For example, `map = {1=>2, 3=>2}` will find all quantum operators with labels 1 and 3, and then map them to 2.

Arguments:

  • g::FeynmanGraph: graph to be modified
  • map: mapping from old labels to the new ones
source
FeynmanDiagram.ComputationalGraphs.remove_all_zero_valued_subgraphs!Method
function remove_all_zero_valued_subgraphs!(graphs::Union{Tuple,AbstractVector{<:AbstractGraph}}; verbose=0)

Recursively removes all zero-valued subgraph(s) in-place in the given graphs.

Arguments:

  • graphs: A collection of graphs to be processed.
  • verbose: Level of verbosity (default: 0).

Returns:

  • Optimized graphs.

source
FeynmanDiagram.ComputationalGraphs.remove_duplicated_leaves!Method
function remove_duplicated_leaves!(graphs::Union{Tuple,AbstractVector{<:AbstractGraph}}; verbose=0, normalize=nothing, kwargs...)

Removes duplicated leaf nodes in-place from a collection of graphs. It also provides optional normalization for these leaves.

Arguments:

  • graphs: A collection of graphs to be processed.
  • verbose: Level of verbosity (default: 0).
  • normalize: Optional function to normalize the graphs (default: nothing).
source
FeynmanDiagram.ComputationalGraphs.remove_zero_valued_subgraphs!Method
function remove_zero_valued_subgraphs!(g::AbstractGraph)

Removes zero-valued (zero subgraph_factor) subgraph(s) of a computational graph `g`. If all subgraphs are zero-valued, the first one (`eldest(g)`) will be retained.

Arguments:

  • g::AbstractGraph: graph to be modified
source
FeynmanDiagram.ComputationalGraphs.remove_zero_valued_subgraphsMethod
function remove_zero_valued_subgraphs(g::AbstractGraph)

Returns a copy of graph `g` with zero-valued (zero subgraph_factor) subgraph(s) removed.
If all subgraphs are zero-valued, the first one (`eldest(g)`) will be retained.

Arguments:

  • g::AbstractGraph: graph to be modified
source
FeynmanDiagram.ComputationalGraphs.replace_subgraph!Method
function replace_subgraph!(g::AbstractGraph, w::AbstractGraph, m::AbstractGraph)

Modifies `g` by replacing the subgraph `w` with a new graph `m`.
For Feynman diagrams, subgraphs `w` and `m` should have the same diagram type, orders, and external indices.

Arguments:

  • g::AbstractGraph: graph to be modified
  • w::AbstractGraph: subgraph to replace
  • m::AbstractGraph: new subgraph
source
FeynmanDiagram.ComputationalGraphs.replace_subgraphMethod
function replace_subgraph(g::AbstractGraph, w::AbstractGraph, m::AbstractGraph)

Creates a modified copy of `g` by replacing the subgraph `w` with a new graph `m`.
For Feynman diagrams, subgraphs `w` and `m` should have the same diagram type, orders, and external indices.

Arguments:

  • g::AbstractGraph: graph to be modified
  • w::AbstractGraph: subgraph to replace
  • m::AbstractGraph: new subgraph
source
FeynmanDiagram.ComputationalGraphs.set_subgraph_factors!Method

function setsubgraphfactors!(g::AbstractGraph, subgraph_factors::AbstractVector, indices::AbstractVector{Int})

Update the specified subgraph factors of graph `g` at indices `indices` to `subgraph_factors`.
By default, calls `set_subgraph_factor!(g, subgraph_factors[i], i)` for each `i` in `indices`.
source
FeynmanDiagram.ComputationalGraphs.set_subgraphs!Method

function set_subgraphs!(g::AbstractGraph, subgraphs::AbstractVector{<:AbstractGraph}, indices::AbstractVector{Int})

Update the specified subgraphs of graph `g` at indices `indices` to `subgraphs`.
By default, calls `set_subgraph!(g, subgraphs[i], i)` for each `i` in `indices`.
source
FeynmanDiagram.ComputationalGraphs.standardize_labels!Method
function standardize_labels!(g::FeynmanGraph)

Finds all labels involved in `g` and its subgraphs and 
modifies `g` by relabeling in standardized order, e.g.,
(1, 4, 5, 7, ...) ↦ (1, 2, 3, 4, ....)

Arguments:

  • g::FeynmanGraph: graph to be relabeled
source
FeynmanDiagram.ComputationalGraphs.standardize_labelsMethod
function standardize_labels!(g::FeynmanGraph)

Finds all labels involved in `g` and its subgraphs and returns 
a copy of `g` relabeled in a standardized order, e.g.,
(1, 4, 5, 7, ...) ↦ (1, 2, 3, 4, ....)

Arguments:

  • g::FeynmanGraph: graph to be relabeled
source
FeynmanDiagram.ComputationalGraphs.subgraph_factorsMethod

function subgraph_factors(g::AbstractGraph, indices::AbstractVector{Int})

Returns the subgraph factors of computational graph `g` at indices `indices`.
By default, calls `subgraph_factor(g, i)` for each `i` in `indices`.
source
FeynmanDiagram.ComputationalGraphs.subgraphsMethod

function subgraphs(g::AbstractGraph, indices::AbstractVector{Int})

Returns the subgraphs of computational graph `g` at indices `indices`.
By default, calls `subgraph(g, i)` for each `i` in `indices`.
source
FeynmanDiagram.ComputationalGraphs.unique_nodes!Function
function unique_nodes!(graphs::AbstractVector{<:AbstractGraph})

Identifies and retrieves unique nodes from a set of graphs.

Arguments:

  • graphs: A collection of graphs to be processed.

Returns:

  • A mapping dictionary from the id of each leaf to the unique leaf node.
source