Reference
Abstract types
SimpleHypergraphs.AbstractHypergraph
— TypeAbstractHypergraph{T} <: AbstractMatrix{T}
An abstract hypergraph type storing information about vertices and hyperedges.
SimpleHypergraphs.AbstractSimpleHypergraph
— TypeAbstractSimpleHypergraph{T} <: AbstractHypergraph{T}
An abstract undirected hypergraph type storing information about vertices and hyperedges.
Creating an undirected hypergraph
SimpleHypergraphs.Hypergraph
— TypeHypergraph{T} <: AbstractSimpleHypergraph{Union{T, Nothing}}
An undirected hypergraph storing information about vertices and hyperedges.
Constructors
Hypergraph{T}(n::Integer,k::Integer) where {T<:Real}
Hypergraph{T,V}(n::Integer, k::Integer;
v_meta=Vector{Union{V,Nothing}}(nothing, n)
) where {T<:Real, V}
Hypergraph{T,E}(n::Integer, k::Integer;
he_meta=Vector{Union{E,Nothing}}(nothing, n)
) where {T<:Real, E}
Hypergraph{T,V,E}(n::Integer, k::Integer;
v_meta=Vector{Union{V,Nothing}}(nothing, n),
he_meta=Vector{Union{E,Nothing}}(nothing, k)
) where {T<:Real, V, E}
Hypergraph{T,V,E,D}(n::Integer, k::Integer,
v_meta=Vector{Union{V,Nothing}}(nothing, n),
he_meta=Vector{Union{E,Nothing}}(nothing, k)
) where {T<:Real,V,E,D<:AbstractDict{Int,T}}
Construct a hypergraph with a given number of vertices and hyperedges. Optionally, values of type V
can be stored at vertices and values of type E
can be stored at hyperedges. By default the hypergraph uses a Dict{Int,T}
for the internal data storage, however a different dictionary such as SortedDict
to ensure result replicability can be used (e.g. when doing stochastic simulations on hypergraphs).
Hypergraph(m::AbstractMatrix{Union{T, Nothing}}) where {T<:Real}
Hypergraph{T,V}(m::AbstractMatrix{Union{T, Nothing}};
v_meta=Vector{Union{V,Nothing}}(nothing, size(m,1))
) where {T<:Real, V}
Hypergraph{T,E}(m::AbstractMatrix{Union{T, Nothing}};
he_meta=Vector{Union{E,Nothing}}(nothing, size(m,2))
) where {T<:Real, E}
Hypergraph{T,V,E}(m::AbstractMatrix{Union{T, Nothing}};
v_meta=Vector{Union{V,Nothing}}(nothing, size(m,1)),
he_meta=Vector{Union{E,Nothing}}(nothing, size(m,2))
) where {T<:Real, V, E}
Hypergraph{T,V,E,D}(m::AbstractMatrix{Union{T, Nothing}};
v_meta=Vector{Union{V,Nothing}}(nothing, size(m,1)),
he_meta=Vector{Union{E,Nothing}}(nothing, size(m,2))
) where {T<:Real, V, E, D<:AbstractDict{Int,T}}
Construct a hypergraph using its matrix representation. In the matrix representation rows are vertices and columns are hyperedges. Optionally, values of type V
can be stored at vertices and values of type E
can be stored at hyperedges. By default the hypergraph uses a Dict{Int,T}
for the internal data storage, however a different dictionary such as SortedDict
to ensure result replicability can be used (e.g. when doing stochastic simulations on hypergraphs).
Hypergraph(g::Graphs.Graph)
Constructs a hypergraph of degree 2 by making a deep copy of Graphs.Graph. A SortedDict
will be used for internal data storage of the hypergraph.
Arguments
T
: type of weight values stored in the hypergraph's adjacency matrixV
: type of values stored in the vertices of the hypergraphE
: type of values stored in the edges of the hypergraphD
: dictionary for storing values the default isDict{Int, T}
n
: number of verticesk
: number of hyperedgesm
: a matrix representation rows are vertices and columns are hyperedgesg
: a graph representation of the hypergraph
SimpleHypergraphs.random_model
— Methodrandom_model(nVertices::Int, nEdges::Int, HType::Type{H}) where {H<:AbstractSimpleHypergraph}
Generate a random undirected hypergraph (in the style of Erdős–Rényi random graphs) without any structural constraints.
The Algorithm
Given two integer parameters nVertices and nEdges (the number of nodes and hyperedges, respectively), the algorithm computes - for each hyperedge he={1,...,m} - a random number s ϵ [1, n] (i.e. the hyperedge size). Then, the algorithm selects uniformly at random s vertices from V to be added in he.
SimpleHypergraphs.random_kuniform_model
— Methodrandom_kuniform_model(nVertices::Int, nEdges::Int, k::Int, HType::Type{H}) where {H<:AbstractSimpleHypergraph}
Generates a k-uniform hypergraph, i.e. an hypergraph where each hyperedge has size k.
The Algorithm
The algorithm proceeds as the random_model, forcing the size of each hyperedge equal to k.
SimpleHypergraphs.random_dregular_model
— Methodrandom_dregular_model(
nVertices::Int,
nEdges::Int,
d::Int,
HType::Type{H}
) where {H<:AbstractSimpleHypergraph}
Generates a d-regular hypergraph, where each node has degree d.
The Algorithm
The algorithm exploits the k-uniform approach described for the randomkuniformmodel method to build a d-regular hypergraph H having nVertices nodes and nEdges edges. It returns the hypergraph H^ dual of H.
SimpleHypergraphs.random_preferential_model
— Methodrandom_preferential_model(
nVertices::Int,
p::Real,
HType::Type{H};
hg::HType = random_model(5,5, HType)
) where {H<:AbstractSimpleHypergraph}
Generate a hypergraph with a preferential attachment rule between nodes, as presented in Avin, C., Lotker, Z., and Peleg, D. Random preferential attachment hyper-graphs. Computer Science 23 (2015).
The Algorithm
The algorithm starts with a random graph with 5 nodes and 5 edges. For this reason, the generated random graph has at least 5 nodes and 5 edges. It iteratively adds a node or a edge, according to a given parameter p, which defines the probability of creating a new node or a new hyperedge.
More in detail, the connections with the new node/hyperedge are generated according to a preferential attachment policy defined by p.
The so-built hypergraph will have nVertices
vertices.
The starting hypergraph can be instantiated as preferred.
Manipulating vertices and hyperedges
SimpleHypergraphs.add_hyperedge!
— Methodadd_hyperedge!(h::Hypergraph{T, V, E, D};
vertices::D = D(), he_meta::Union{E,Nothing}=nothing
) where {T <: Real, V, E, D <: AbstractDict{Int,T}}
Adds a hyperedge to a given undirected hypergraph h
. Optionally, existing vertices can be added to the created hyperedge. The paramater vertices
represents a dictionary of vertex identifiers and values stored at the hyperedges. Additionally, a value can be stored with the hyperedge using the he_meta
keyword parameter.
SimpleHypergraphs.add_vertex!
— Methodadd_vertex!(h::Hypergraph{T, V, E, D};
hyperedges::D = D(), v_meta::Union{V,Nothing} = nothing
) where {T <: Real, V, E, D <: AbstractDict{Int,T}}
Adds a vertex to a given undirected hypergraph h
. Optionally, the vertex can be added to existing hyperedges. The hyperedges
parameter presents a dictionary of hyperedge identifiers and values stored at the hyperedges. Additionally, a value can be stored with the vertex using the v_meta
keyword parameter.
SimpleHypergraphs.set_vertex_meta!
— Methodset_vertex_meta!(h::Hypergraph{T, V, E, D}, new_value::Union{V,Nothing},
id::Int) where {T <: Real, V, E, D <: AbstractDict{Int,T}}
Sets a new meta value new_value
for the vertex id
in the hypergraph h
.
SimpleHypergraphs.get_vertex_meta
— Methodget_vertex_meta(h::Union{Hypergraph{T, V, E, D}, DirectedHypergraph{T, V, E, D}}, id::Int
) where {T <: Real, V, E, D <: AbstractDict{Int,T}}
Returns a meta value stored at the vertex id
in the undirected hypergraph h
.
SimpleHypergraphs.set_hyperedge_meta!
— Methodset_hyperedge_meta!(h::Hypergraph{T, V, E, D},
new_value::Union{E,Nothing}, id::Int
) where {T <: Real, V, E, D <: AbstractDict{Int,T}}
Sets a new meta value new_value
for the hyperedge id
in the undirected hypergraph h
.
SimpleHypergraphs.get_hyperedge_meta
— Methodget_hyperedge_meta(h::Hypergraph{T, V, E, D}, id::Int)
where {T <: Real, V, E, D <: AbstractDict{Int,T}}
Returns a meta value stored at the hyperedge id
in the undirected hypergraph h
.
SimpleHypergraphs.remove_vertex!
— Methodremove_vertex!(h::Hypergraph, v::Int)
Removes the vertex v
from a given undirected hypergraph h
. Note that running this function will cause reordering of vertices in the hypergraph; the vertex v
will replaced by the last vertex of the hypergraph and the list of vertices will be shrunk.
SimpleHypergraphs.remove_hyperedge!
— Methodremove_hyperedge!(h::Hypergraph, e::Int)
Removes the hyperedge e
from a given undirected hypergraph h
. Note that running this function will cause reordering of hyperedges in the hypergraph: the hyperedge e
will replaced by the last hyperedge of the hypergraph and the list of hyperedges (and hyperedge metadata) will be shrunk.
Hypergraph array getters and setters
Normally you work with a hypergraph via array setters, for example the code below craete an Hypergraph and add vertex one to hyperedges 2 and 3 with weight 5:
h = Hypergraph{Int64}(2,3);
h[1, 2:3] .= 5;
h
# output
2×3 Hypergraph{Int64, Nothing, Nothing, Dict{Int64, Int64}}:
nothing 5 5
nothing nothing nothing
Base.getindex
— MethodBase.getindex(h::H, idx::Vararg{Int,2}) where {H <: AbstractSimpleHypergraph}
Returns a value for a given vertex-hyperedge pair idx
for an undirected hypergraph h
. If a vertex does not belong to a hyperedge nothing
is returned.
Base.setindex!
— MethodBase.setindex!(h::H, ::Nothing, idx::Vararg{Int,2}) where {H <: AbstractSimpleHypergraph}
Removes a vertex from a given hyperedge for an undirected hypergraph h
and a given vertex-hyperedge pair idx
. Note that trying to remove a vertex from a hyperedge when it is not present will not throw an error.
Base.setindex!
— MethodBase.setindex!(h::H, v::Real, idx::Vararg{Int,2}) where {H <: AbstractSimpleHypergraph}
Adds a vertex to an undirected hyperedge (represented by indices idx
) and assigns value v
to be stored with that assignment.
Hypergraph representation as Graphs.jl simple graphs
The goal of those methods is to provide a way to manipulate a hypergraph using the methods from the Graphs.jl library. This has been achieved by providing types that are subtypes of the Graphs.SimpleGraphs.AbstractSimpleGraph{Int}
type along with appropiate methods.
SimpleHypergraphs.BipartiteView
— TypeBipartiteView{T<:Real} <: Graphs.SimpleGraphs.AbstractSimpleGraph{Int}
Represents a bipartite view of a hypergraph h
. Note this is a view - changes to the original hypergraph will be automatically reflected in the view.
Constructors
BipartiteView(::H) where {H<:AbstractHypergraph}
The bipartite view of a hypergraph is suitable for processing with the Graphs.jl package. Several Graphs methods are provided for the compability.
SimpleHypergraphs.shortest_path
— Methodshortest_path(b::BipartiteView{H}, source::Int, target::Int) where {H<:AbstractSimpleHypergraph}
Finds a single shortest path in a graph b
between vertices source
and target
. Note that if several paths of the same length exist, only one will be returned.
SimpleHypergraphs.TwoSectionView
— TypeTwoSectionView{H<:AbstractHypergraph} <: Graphs.SimpleGraphs.AbstractSimpleGraph{Int64}
Represents a 2-section view of a hypergraph h
. Note (1) this is a view - changes to the original hypergraph will be automatically reflected in the view.
Note (2) The view will only work correctly for hypergraphs not having overlapping hyperedges. To check whether a graph has overlapping edges try has_overlapping_hedges(h)
- for such graph you need to fully materialize it rather than use a view. This can be achieved via the get_twosection_adjacency_mx(h)
method.
Constructors
TwoSectionView(::H<:AbstractHypergraph)
The 2-section view of a hypergraph is suitable for processing with the Graphs.jl package. Several Graphs methods are provided for the compability.
SimpleHypergraphs.shortest_path
— Methodshortest_path(t::TwoSectionView{H}, source::Int, target::Int) where {H<:AbstractSimpleHypergraph}
Finds a single shortest path in a graph b
between vertices source
and target
. Note that if several paths of the same length exist, only one will be returned.
Hypergraph info
Base.size
— MethodBase.size(h::AbstractHypergraph)
Returns the size of hypergraph h
. The result is a tuple of the number of vertices and the number of hyperedges
SimpleHypergraphs.nhv
— Methodnhv(h::H) where {H <: AbstractSimpleHypergraph}
Return the number of vertices in the undirected hypergraph h
.
SimpleHypergraphs.nhe
— Methodnhe(h::H) where {H <: AbstractSimpleHypergraph}
Return the number of hyperedges in the undirected hypergraph h
.
SimpleHypergraphs.getvertices
— Methodgetvertices(h::H, he_id::Int) where {H <: AbstractSimpleHypergraph}
Returns vertices from an undirected hypergraph a
for a given hyperedge he_id
.
SimpleHypergraphs.gethyperedges
— Methodgethyperedges(h::H, v_id::Int) where {H <: AbstractSimpleHypergraph}
Returns hyperedges for a given vertex v_id
in an undirected hypergraph h
.
SimpleHypergraphs.get_connected_components
— Methodget_connected_components(h::H) where {H <: AbstractSimpleHypergraph}
Return an array of connected components in the undirected hypergraph h
(array of vectors of vertices) using recurrence.
SimpleHypergraphs.conductance
— Methodconductance(h::H, subset::Set{Int})::Float64 where {H<:AbstractSimpleHypergraph}
Calculate unweighted hypergraph conductance of subset
. note: ∅ ⊊ subset
⊊ 1:nhv(h)
For more information see 1. Introduction
at: Generalizing the Hypergraph Laplacian via a Diffusion Process with Mediators, authors: T-H. Hubert Chan, Xhibin Liang.
SimpleHypergraphs.get_twosection_adjacency_mx
— Methodget_twosection_adjacency_mx(h::H; count_self_loops::Bool=false,
replace_weights::Union{Nothing,Real}=nothing ) where {H<:AbstractSimpleHypergraph}
Returns an adjacency matrix for a two section view of a hypergraph h
.
SimpleHypergraphs.random_walk
— Methodrandom_walk(
h::H,
start::Int;
heselect::Function,
vselect::Function,
) where {H <: AbstractSimpleHypergraph}
Return a next vertex visited in assuming a random walk starting from vertex start
. First a hyperedge is sampled with weights proportional to heselect
function (by default each hyperedge is sampled with the same probability). Next a vertex within hyperedge is with weights proportional to vselect
function (by default each vertex, including the source, is sampled with the same probability).
heselect
and vselect
functions take two arguments a Hypergraph
and respectively a vertex identifier or a hyperedge identifier. The return values of both functions should be respectively a list of hyperedges or vertices and their weights.
SimpleHypergraphs.dual
— Methoddual(h::Hypergraph)
Return the dual of the hypergraph h
.
NOTE h
needs to have at least one dimension greater than 0.
Graphs.modularity
— MethodGraphs.modularity(h::H, partition::Vector{Set{Int}},
ha::HypergraphAggs=HypergraphAggs(h)) where {H<:AbstractSimpleHypergraph}
Calculates the strict modularity of a hypergraph h
for a given partition
using the precomputed aggregates ha
.
SimpleHypergraphs.HypergraphAggs
— TypeHypergraphAggs(h::H) where {H<:AbstractSimpleHypergraph}
Precomputes vertex and edge basic stats for a hypergraph. The stats are being used for efficiency reasons by community search algorithms.
SimpleHypergraphs.randompartition
— Methodrandompartition(N::Int, n::Int)::Vector{Set{Int}}
Generates a random partition for graph having N
vertices into n
subsets.
SimpleHypergraphs.randompartition
— Methodrandompartition(h::Hypergraph, n::Int)::Vector{Set{Int}}
Generates a random partition for vertices of a hypergraph h
into n
subsets.
SimpleHypergraphs.AbstractCommunityFinder
— TypeThe base type for all algorithms representing various community search patterns.
SimpleHypergraphs.CFModularityRandom
— TypeCFModularityRandom(n::Int, reps::Int) <: AbstractCommunityFinder
Represents a random search over the hypergraph h
that finds a partition into n
communities (subsets) having the maximum modularity value. During the search reps
random n
-partitions will be evaluated. If there are many partitions having the same value the first one that was randomly found will be returned.
SimpleHypergraphs.CFModularityCNMLike
— TypeCFModularityCNMLike(n::Int, reps::Int) <: AbstractCommunityFinder
Represents a CNM-Like algorithm for finding communities. In the algorithm we start with a partition where each node is in its own part. Then in each step, we randomly select a hyperedge. Subsequently, we consider merging each set of that parts it touches. We actually merge the parts if the new best modularity is at least as high as the modularity from the previous step. The algortithm iterates through reps
of repetitions.
For more information see Algorithm 1
at: Clustering via Hypergraph Modularity (submitted to Plos ONE), auhtors: Bogumił Kamiński, Valerie Poulin, Paweł Prałat, Przemysław Szufel, Francois Theberge
SimpleHypergraphs.findcommunities
— Methodfindcommunities(h::H, method::CFModularityRandom) where {H<:AbstractSimpleHypergraph}
Makes a random search over the hypergraph h
and finds a partition into method.n
communities (subsets) having the maximum modularity value. During the search method.reps
random n
-partitions will be evaluated. If there are many partitions having the same value the first one that was randomly found will be returned.
Returns a NamedTuple
where the field bp
contains partition and the field bm
contains the modularity value for that partition.
SimpleHypergraphs.findcommunities
— Methodfindcommunities(h::H, method::CFModularityCNMLike) where {H<:AbstractSimpleHypergraph}
Iterates a CNM-Like algorithm for finding communities. In the algorithm we start with a partition where each node is in its own part. Then in each step, we randomly select a hyperedge. Subsequently, we consider merging each set of that parts it touches. We actually merge the parts if the new best modularity is at least as high as the modularity from the previous step.
Returns a NamedTuple
where the field bp
contains partition and the field bm
contains the modularity value for that partition, finally, the fiel mod_history
represents modularities achieved in subsequent steps of the algorithm.
For more information see Algorithm 1
at: Clustering via Hypergraph Modularity (submitted to Plos ONE), authors: Bogumił Kamiński, Valerie Poulin, Paweł Prałat, Przemysław Szufel, Francois Theberge
I/O
Undirected hypergraphs can be saved as and loaded from JSON- and HGF-formatted files.
SimpleHypergraphs.hg_save
— Function`` hgsave(io::IO, h::H, format::HGFFormat) where {H <: AbstractSimpleHypergraph}
Saves an undirected hypergraph h
to an output stream io
in hgf
format.
TODO: what to do about metadata?
hg_save(io::IO, h::Hypergraph, format::JSON_Format)
Saves an undirected hypergraph h
to an output stream io
in json
format.
If h
has Composite Types
either for vertex metadata or hyperedges metadata, the user has to explicit tell the JSON3 package about it, for instance using:
JSON3.StructType(::Type{MyType}) = JSON3.Struct()
.
See the (JSON3.jl documentation)[https://github.com/quinnj/JSON3.jl] for more details.
The json
in output contains the following information (keys):
n
: number of verticesk
: number of hyperedgesm
: a matrix representation ofh
where rows are vertices and columns are hyperedgesv2he
: mapping vertices to hyperedgesv_meta
: vertices metadatahe_meta
: hyperedges metadata
hg_save(
fname::AbstractString, h::AbstractHypergraph;
format::Abstract_HG_format=HGF_Format()
)
Saves a hypergraph h
to a file fname
in the specified format
. The default saving format is hgf
.
SimpleHypergraphs.hg_load
— Functionhg_load(
io::IO,
format::HGF_Format;
HType::Type{H} = Hypergraph,
T::Type{U} = Bool,
D::Type{<:AbstractDict{Int, U}} = Dict{Int, T},
) where {U <: Real, H <: AbstractSimpleHypergraph}
Loads a hypergraph from a stream io
from hgf
format.
Arguments
T
: type of weight values stored in the hypergraph's adjacency matrixD
: dictionary for storing values the default isDict{Int, T}
Skips a single initial comment.
hg_load(
io::IO,
T::Type{H},
format::JSON_Format;
T::Type{U} = Bool,
D::Type{<:AbstractDict{Int, U}} = Dict{Int,U},
V = Nothing,
E = Nothing
) where {H <: AbstractHypergraph, U <: Real}
Loads a hypergraph from a stream io
from json
format.
Arguments
T
: type of weight values stored in the hypergraph's adjacency matrixD
: dictionary for storing values the default isDict{Int, T}
V
: type of values stored in the vertices of the hypergraphE
: type of values stored in the edges of the hypergraph
hg_load(
fname::AbstractString;
format::Abstract_HG_format = HGF_Format(),
HType::Type{H} = Hypergraph,
T::Type{U} = Bool,
D::Type{<:AbstractDict{Int, U}} = Dict{Int, T},
V = Nothing,
E = Nothing
) where {U <: Real, H <: AbstractSimpleHypergraph}
Loads a hypergraph from a file fname
. The default saving format is hgf
.
Arguments
HType
: type of hypergraph to store data inT
: type of weight values stored in the hypergraph's adjacency matrixD
: dictionary for storing values the default isDict{Int, T}
V
: type of values stored in the vertices of the hypergraphE
: type of values stored in the edges of the hypergraph
Hypergraph Visualization
Currently, visualization is only supported for undirected hypergraphs.
SimpleHypergraphs.draw
— Functionfunction draw(
h::H,
type::Type{GraphBased};
element::Union{String, Int}=get_next_div_id(),
width::Int=500,
height::Int=500,
radius::Real=10,
node_radii::Union{AbstractVector{<:Real}, Nothing}=nothing,
node_color::String="#999",
node_colors::Union{AbstractVector{String}, Nothing}=nothing,
node_stroke::Union{String, Nothing} = nothing,
node_strokes::Union{AbstractVector{String}, Nothing}=nothing,
stroke_width::Real=0,
stroke_widths::Union{AbstractVector{<:Real}, Nothing}=nothing,
node_opacity::Real=1,
node_opacities::Union{AbstractVector{<:Real}, Nothing}=nothing,
stroke_opacity::Real=1,
stroke_opacities::Union{AbstractVector{<:Real}, Nothing}=nothing,
with_node_labels::Bool=false,
node_labels::Union{AbstractVector{String}, Nothing}=nothing,
with_node_metadata_hover::Bool=false,
with_node_weight::Bool=false,
he_colors::Union{AbstractVector{String}, Nothing}=nothing,
with_he_labels::Bool=false,
he_labels::Union{AbstractVector{String}, Nothing}=nothing,
with_he_metadata_hover::Bool=false
) where {H<:AbstractSimpleHypergraph}
Draw a hypergraph h
in a web-based environment (e.g. Jupyter Notebook), using a js script based on the library (D3)[https://d3js.org/]. Each hyperedge he
is represented by a fake vertex fv
to which each vertex v ∈ he
is connected.
Arguments
h
: the hypergraph to drawtype
: how the hypergraph will be drawn. Iftype=GraphBased
, each hyperedge
will be represented as a vertex (see above)
width
: width of the figureheight
: height of the figureradius
: same default radius for each vertex (represented as a circle)node_radii
: distinct radius values for each vertexnode_color
: same default color for each vertexnode_colors
: distinct node colors for each vertexnode_stroke
: same default stroke for each vertexnode_strokes
: distinct node strokes for each vertexstroke_width
: same default stroke-width for each vertexstroke_widths
: distinct stroke-width values for each vertexnode_opacity
: same default opacity for each vertexnode_opacities
: distinct node-opacity values for each vertexstroke_opacity
: same default stroke-opacity for each vertexstroke_opacities
: distinct stroke-opacity values for each vertexwith_node_labels
: whether displaying node labelsnode_labels
: node labels to be shownwith_node_metadata_hover
: whether displaying node metadata when hovering each vertexwith_node_weight
: whether displaying node weights within each hyperedgehe_colors
: distinct hyperedge colors for each hyperedgewith_he_labels
: whether displaying hyoeredges labelswith_he_metadata_hover
: whether displaying hyperedge metadata when hovering each hyperedge
draw(
h::H,
type::Type{HyperNetX};
width::Int=10,
height::Int=10,
node_labels::Union{Dict{Int, String}, Nothing}=nothing,
edge_labels::Union{Dict{Int, String}, Nothing}=nothing,
collapse_nodes::Bool=false,
collapse_edges::Bool=false,
pos::Union{Dict{Int,Pair{Int,Int}}, Nothing}=nothing,
with_color::Bool=true,
with_node_counts::Bool=false,
with_edge_counts::Bool=false,
layout::PyObject=nx.spring_layout,
layout_kwargs::Dict=Dict{String, Any}(),
ax::Union{PyObject, Nothing}=nothing,
no_border::Bool=false,
edges_kwargs::Dict=Dict{String, Any}(),
nodes_kwargs::Dict=Dict{String, Any}(),
edge_labels_kwargs::Dict=Dict{String, Any}(),
node_labels_kwargs::Dict=Dict{String, Any}(),
with_edge_labels::Bool=true,
with_node_labels::Bool=true,
label_alpha::Float64=.35
) where {H<:AbstractSimpleHypergraph}
Draw a hypergraph h
as an Euler diagram, using the library HyperNetX.
Arguments
h
: the hypergraph to drawtype
: how the hypergraph will be drawn. Iftype=HyperNetX
, the hypergraph will be represented as a Euler Diagramwidth
: width of the figureheight
: height of the figurenode_labels
: node labels to be shownedge_labels
: edge labels to be showncollapse_nodes
: draws the hypergraph gotten by identifying nodes contained by the same edges (from HyperNetX)collapse_edges
: draws the hypergraph gotten by identifying edges containing the same nodes (from HyperNetX)no_border
: indicates wheter the figure should have a border
For more details about the other parameters, please refer to the library HyperNetX.