API Documentation
OptiGraph
Plasmo.OptiGraph
— TypeOptiGraph()
Create an empty OptiGraph. An OptiGraph extends a JuMP.AbstractModel and supports most JuMP.Model functions.
Plasmo.@optinode
— Macro@optinode(optigraph, expr...)
Add a new optinode to optigraph
. The expression expr
can either be
- of the form
nodename
creating a single optinode with the variable namevarname
- of the form
nodename[...]
or[...]
creating a container of optinodes using JuMP Containers
Plasmo.@linkconstraint
— Macro@linkconstraint(graph::OptiGraph, expr)
Add a linking constraint described by the expression expr
.
@linkconstraint(graph::OptiGraph, ref[i=..., j=..., ...], expr)
Add a group of linking constraints described by the expression expr
parametrized by i
, j
, ...
The @linkconstraint macro works the same way as the JuMP.@constraint
macro.
Plasmo.optigraph_reference
— Functionoptigraph_reference(graph::OptiGraph)::OptiGraph
Create a new optigraph with the same optinodes and optiedges as graph
. Useful for defining an optigraph over existing nodes and edges without recreating them. Note that any changes to the optinodes and optiedges in the returned optigraph will take place in the original graph
.
Plasmo.add_subgraph!
— Functionadd_subgraph!(graph::OptiGraph, subgraph::OptiGraph)::OptiGraph
Add the sub-optigraph subgraph
to the higher level optigraph graph
. Returns the original graph
Plasmo.subgraph
— Functionsubgraph(graph::OptiGraph, idx::Int64)
Retrieve the the subgraph in graph
at index idx
.
Plasmo.subgraphs
— Functionsubgraphs(optigraph::OptiGraph)::Vector{OptiGraph}
Retrieve the local subgraphs of optigraph
.
Plasmo.all_subgraphs
— Functionall_subgraphs(graph::OptiGraph)::Vector{OptiGraph}
Retrieve all of the contained subgraphs of graph
, including nested subgraphs. The order of the subgraphs in the returned vector starts with the local subgraphs in optigraph
and then appends the nested subgraphs for each local subgraph.
Plasmo.subgraph_by_index
— Functionsubgraph_by_index(graph::OptiGraph, index::Int64)::OptiGraph
Recursively search optigraph graph
for the subngraph at index
by traversing subgraphs. Note that the subgraph is not unique to the index. Since the search is depth-first, the subgraph returned may be different if the overall graph
structure changes.
Plasmo.num_subgraphs
— Functionnum_subgraphs(graph::OptiGraph)::Int64
Retrieve the number of local subgraphs in graph
.
Plasmo.num_all_subgraphs
— Functionnum_all_subgraphs(graph::OptiGraph)::Int64
Retrieve the number of subgraphs in graph
including nested subgraphs.
Plasmo.has_subgraphs
— Functionhas_subgraphs(graph::OptiGraph)::Bool
Check whether graph
contains subgraphs.
Plasmo.add_node!
— Functionadd_node!(graph::OptiGraph)::OptiNode
Create a new OptiNode
and add it to graph
. Returns the added optinode.
add_node!(graph::OptiGraph, m::JuMP.Model)::OptiNode
Add a new optinode to graph
and set its model to the JuMP.Model
m
.
add_node!(graph::OptiGraph, optinode::OptiNode)::OptiNode
Add the existing optinode
(Created with OptiNode()
) to the optigraph graph
.
Plasmo.optinode
— Functionoptinode(m::JuMP.Model)
Retrieve the optinode corresponding to JuMP model m
optinode(var::JuMP.VariableRef)
Retrieve the optinode corresponding to JuMP VariableRef
optinode(var::JuMP.VariableRef)::OptiNode
Retrieve the optinode corresponding to JuMP ConstraintRef
optinode(graph::OptiGraph, index::Int64)
Retrieve the local optinode in graph
at index
. This does not retrieve nodes in subgraphs.
Plasmo.optinodes
— Functionoptinodes(graph::OptiGraph)::Vector{OptiNode}
Retrieve the optinodes in graph
. Note that this returns the local optinodes contained directly in graph
and excludes nodes contained in subgraphs of graph
.
Plasmo.all_nodes
— Functionall_nodes(graph::OptiGraph)::Vector{OptiNode}
Recursively collect all optinodes in graph
by traversing each of its subgraphs.
Plasmo.optinode_by_index
— Functionoptinode_by_index(graph::OptiGraph, index::Int64)::OptiNode
Recursively search optigraph graph
for the optinode at index
by traversing subgraphs. Note that the optinode is not unique to the index. Since the search is depth-first, the optinode returned may be different if the subgraph structure changes.
Plasmo.add_optiedge!
— Functionadd_optiedge!(graph::OptiGraph, optinodes::Vector{OptiNode})::OptiEdge
Add an optiedge to optigraph graph
that connects optinodes
. If edge already exists, return it.
Plasmo.optiedge
— Functionoptiedge(graph::OptiGraph, index::Int64)
Retrieve the local optiedge in graph
at index
optiedge(graph::OptiGraph, nodes::OrderedSet{OptiNode})
Retrieve the optiedge in graph
that connects the optinodes in the OrderedSet of nodes
.
optiedge(graph::OptiGraph, nodes::OptiNode...)
Retrieve the optiedge in graph
that connects nodes
.
Plasmo.optiedges
— Functionoptiedges(graph::OptiGraph)::Vector{OptiEdge}
Retrieve the local optiedges in graph
.
Plasmo.all_edges
— Functionall_edges(graph::OptiGraph)::Vector{OptiEdge}
Recursively collect all optiedges in graph
by traversing each of its subgraphs.
Plasmo.optiedge_by_index
— Functionoptiedge_by_index(graph::OptiGraph, index::Int64)::OptiEdge
Recursively search optigraph graph
for the edge at index
by traversing subgraphs. Note that the edge is not unique to the index. Since the search is depth-first, the optiedge returned may be different if the subgraph structure changes.
Plasmo.linkconstraints
— Functionlinkconstraints(graph::OptiGraph)::Vector{LinkConstraintRef}
Retrieve the local linking constraints in graph
. Returns a vector of the linking constraints.
Plasmo.all_linkconstraints
— Functionall_linkconstraints(graph::OptiGraph)::Vector{LinkConstraintRef}
Recursively collect all linkconstraints in graph
by traversing each of its subgraphs.
Plasmo.num_linkconstraints
— Functionnum_linkconstraints(node::OptiNode)
Return the number of link-constraints incident to optinode node
num_linkconstraints(graph::OptiGraph)::Int64
Retrieve the number of local linking constraints in graph
. Does not include linkconstraints in subgraphs.
Plasmo.has_objective
— Functionhas_objective(node::OptiNode)
Check whether optinode node
has a non-empty linear or quadratic objective function
has_objective(graph::OptiGraph)::Bool
Check whether optigraph graph
has an affine or quadratic objective function set.
Plasmo.has_nl_objective
— Functionhas_nl_objective(node::OptiNode)
Check whether optinode node
has a nonlinear objective function
has_nl_objective(graph::OptiGraph)::Bool
Check whether any optinode in graph
has a nonlinear objective function.
Base.getindex
— MethodBase.getindex(graph::OptiGraph, node::OptiNode)
Retrieve the index of the optinode node
in graph
.
Base.getindex
— MethodBase.getindex(graph::OptiGraph, optiedge::OptiEdge)::Int64
Retrieve the index of the optiedge
in graph
.
OptiNode
Plasmo.OptiNode
— TypeOptiNode()
Creates an empty OptiNode
. Does not add it to an OptiGraph
.
Plasmo.jump_model
— Functionjump_model(node::OptiNode)
Get the underlying JuMP.Model
from the optinode node
.
Plasmo.set_model
— Functionset_model(node::OptiNode, m::AbstractModel)
Set the JuMP model m
to optinode node
.
Plasmo.label
— Functionlabel(node::OptiNode)
Get the label for optinode node
.
Plasmo.set_label
— Functionsetlabel(node::OptiNode, label::Symbol)
Set the label for optinode node
to label
. This is what gets printed.
Plasmo.is_node_variable
— Functionis_node_variable(node::OptiNode, vref::JuMP.AbstractVariableRef)
Checks whether the variable reference vref
belongs to the optinode node
.
is_node_variable(vref::JuMP.AbstractVariableRef)
Checks whether the variable reference vref
belongs to any OptiNode
.
Plasmo.is_set_to_node
— Functionis_set_to_node(m::JuMP.AbstractModel)
Checks whether the JuMP model m
is set to any OptiNode
Base.getindex
— MethodBase.getindex(node::OptiNode, symbol::Symbol)
Support retrieving node attributes via symbol lookup. (e.g. node[:x])
Base.setindex
— MethodBase.setindex(node::OptiNode, value::Any, symbol::Symbol)
Support retrieving node attributes via symbol lookup. (e.g. node[:x])
Plasmo.num_linked_variables
— Functionnum_linked_variables(node::OptiNode)
Return the number of node variables on node
that are linked to other nodes
OptiEdge
Plasmo.OptiEdge
— TypeOptiEdge
The OptiEdge
type. Typically created from @linkconstraint
. Contains the set of its supporting optionodes, as well as references to its underlying linking constraints.
Plasmo.LinkConstraint
— TypeLinkConstraint{F <: JuMP.AbstractJuMPScalar,S <: MOI.AbstractScalarSet} <: AbstractLinkConstraint
Type inherits JuMP.AbstractConstraint. Contains a func and set used to describe coupling between optinodes.
LinkConstraint(con::JuMP.ScalarConstraint)
Creates a linking constraint from a JuMP.ScalarConstraint.
LinkConstraint(ref::LinkConstraintRef)
Retrieves a linking constraint from a LinkConstraintRef.
Plasmo.LinkConstraintRef
— TypeLinkConstraintRef
A constraint reference to a linkconstraint. Stores linkconstraint id and the optiedge it belong to.
Plasmo.attached_node
— Functionattached_node(con::LinkConstraint)
Retrieve the attached node on linkconstraint con
Plasmo.set_attached_node
— Functionset_attached_node(con::LinkConstraint,node::OptiNode).
Set the linkconstraint con
to optinode node
. Mostly useful for algorithms that need an "owning" node on a linkconstraint
Extended Functions
JuMP.all_variables
— FunctionJuMP.all_variables(node::OptiNode)::Vector{JuMP.VariableRef}
Retrieve all of the variables on the optinode node
.
JuMP.all_variables(graph::OptiGraph)::Vector{JuMP.VariableRef}
Retrieve a list of all variables in optigraph graph.
JuMP.set_optimizer
— FunctionJuMP.set_optimizer(graph::OptiGraph, optimizer_constructor::Any)
Set an MOI optimizer onto the optigraph graph
. Works exactly the same as using JuMP to set an optimizer.
Example
graph = OptiGraph()
set_optimizer(graph, GLPK.Optimizer)
JuMP.optimize!
— FunctionJuMP.optimize!(graph::OptiGraph)
Optimize the optigraph graph
with the current set optimizer
Example
graph = OptiGraph()
set_optimizer(graph, GLPK.Optimizer)
optimize!(graph)
JuMP.objective_function
— FunctionJuMP.objective_function(node::OptiNode)
Retrieve the objective function on optinode node
JuMP.objective_function(graph::OptiGraph)
Retrieve the current graph objective function.
JuMP.value
— FunctionJuMP.value(node::OptiNode, vref::VariableRef)
Get the variable value of vref
on the optinode node
. This value is always the local node value, not the value the node variable takes when solved as part of a larger OptiGraph
.
JuMP.value(graph::OptiGraph, vref::VariableRef)
Get the variable value of vref
on the optigraph graph
. This value corresponds to the optinode variable value obtained by solving graph
which contains said optinode.
JuMP.dual
— FunctionJuMP.dual(c::JuMP.ConstraintRef{OptiNode,NonlinearConstraintIndex})
Get the dual value on a nonlinear constraint on an OptiNode
JuMP.dual(graph::OptiGraph, linkref::LinkConstraintRef)
Retrieve the dual value of linkref
on optigraph graph
.
JuMP.num_variables
— FunctionJuMP.num_variables(node::OptiNode)
Get the number of variables on optinode node
JuMP.num_variables(graph::OptiGraph)::Int64
Retrieve the number of local node variables in graph
. Does not include variables in subgraphs.
JuMP.num_constraints
— FunctionJuMP.num_constraints(node::OptiNode)
Get the number of constraints on optinode node
JuMP.num_constraints(graph::OptiGraph)
Retrieve the number of local node constraints in graph
. Does not include constraints in subgraphs.
JuMP.object_dictionary
— FunctionJuMP.object_dictionary(node::OptiNode)
Get the underlying object dictionary of optinode node
JuMP.object_dictionary(graph::OptiGraph)
Retrieve the object dictionary of optigraph graph
JuMP.add_variable
— FunctionJuMP.add_variable(node::OptiNode, v::JuMP.AbstractVariable, name::String="")
Add variable v
to optinode node
. This function supports use of the @variable
JuMP macro. Optionally add a base_name
to the variable for printing.
JuMP.add_constraint
— FunctionJuMP.add_constraint(node::OptiNode, con::JuMP.AbstractConstraint, base_name::String="")
Add a constraint con
to optinode node
. This function supports use of the @constraint JuMP macro.
JuMP.add_nonlinear_constraint
— FunctionJuMP.add_nonlinear_constraint(node::OptiNode, expr::Expr)
Add a non-linear constraint to an optinode using a Julia expression.
JuMP.num_nonlinear_constraints
— FunctionJuMP.num_nonlinear_constraints(node::OptiNode)
Get the number of nonlinear constraints on optinode node
JuMP.list_of_constraint_types
— FunctionJuMP.list_of_constraint_types(node::OptiNode)
Get a list of constraint types on optinode node
JuMP.list_of_constraint_types(graph::OptiGraph)
Retrieve a list of the constraint types in optigraph graph
JuMP.all_constraints
— FunctionJuMP.all_constraints(node::OptiNode,F::DataType,S::DataType)
Get all constraints on optinode node
of function type F
and set S
JuMP.all_constraints(graph::OptiGraph, F::DataType, S::DataType)
Retrieve a list of contraints with function type F
and set S
in optigraph graph
JuMP.objective_value
— FunctionJuMP.objective_value(graph::OptiGraph)
Retrieve the current objective value on optigraph graph
.
JuMP.objective_sense
— FunctionJuMP.objective_function(graph::OptiGraph)::MOI.OptimizationSense
Retrieve the current graph objective sense.
JuMP.set_objective
— FunctionJuMP.set_objective(graph::OptiGraph, sense::MOI.OptimizationSense, func::JuMP.AbstractJuMPScalar)
Set the objective of graph
to the optimization sense
and func
.
JuMP.set_nonlinear_objective
— FunctionJuMP.set_nonlinear_objective(optinode::OptiNode, sense::MOI.OptimizationSense, obj::Any)
Set a nonlinear objective on optinode node
JuMP.set_objective_function
— FunctionJuMP.set_objective_function(graph::OptiGraph, x::JuMP.VariableRef)
Set a single variable objective function on optigraph graph
JuMP.set_objective_function(graph::OptiGraph, expr::JuMP.GenericAffExpr)
Set an affine objective function on optigraph graph
JuMP.set_objective_function(graph::OptiGraph, expr::JuMP.GenericQuadExpr)
Set a quadratic objective function on optigraph graph
JuMP.set_objective_sense
— FunctionJuMP.set_objective_sense(graph::OptiGraph, sense::MOI.OptimizationSense)
Set the current graph objective sense to sense
.
JuMP.NLPEvaluator
— FunctionJuMP.NLPEvaluator(node::OptiNode)
Retrieve the underlying JuMP NLP evaluator on optinode node
JuMP.termination_status
— FunctionJuMP.termination_status(node::OptiNode)
Return the termination status on optinode node
JuMP.termination_status(graph::OptiGraph)
Retrieve the current termination status of optigraph graph
.
Graph Processing and Partitioning
Plasmo.graph_backend
— Functiongraph_backend(graph::OptiGraph)
Retrieve the underlying hypergraph backend of an optigraph.
Plasmo.HyperGraph
— TypeHyperGraph
A simple hypergraph type. Contains attributes for vertices and hyperedges.
Plasmo.hyper_graph
— Functionhyper_graph(graph::OptiGraph)
Retrieve a hypergraph representation of the optigraph graph
. Returns a HyperGraph
object, as well as a dictionary that maps hypernodes and hyperedges to the original optinodes and optiedges.
LightGraphs.SimpleGraphs.clique_graph
— Functionclique_graph(graph::OptiGraph)
Retrieve a standard graph representation of the optigraph graph
. Returns a LightGraphs.Graph
object, as well as a dictionary that maps vertices and edges to the optinodes and optiedges.
Plasmo.edge_graph
— Functionedge_graph(optigraph::OptiGraph)
Retrieve the edge-graph representation of optigraph
. This is sometimes called the line graph of a hypergraph. Returns a LightGraphs.Graph
object, as well as a dictionary that maps vertices and edges to the optinodes and optiedges.
Plasmo.edge_hyper_graph
— Functionedge_hyper_graph(graph::OptiGraph)
Retrieve an edge-hypergraph representation of the optigraph graph
. Returns a HyperGraph
object, as well as a dictionary that maps hypernodes and hyperedges to the original optinodes and optiedges. This is also called the dual-hypergraph representation of a hypergraph.
Plasmo.bipartite_graph
— Functionbipartite_graph(optigraph::OptiGraph)
Create a bipartite graph representation from optigraph
. The bipartite graph contains two sets of vertices corresponding to optinodes and optiedges respectively.
Plasmo.Partition
— TypePartition(hypergraph::HyperGraph,node_membership_vector::Vector{Int64},ref_map::Dict)
Create a partition of optinodes using hypergraph
, node_membership_vector
, and 'refmap'. The 'refmap' is a dictionary that maps hypernode indices (integers) and hyperedge indices (tuples) back to optinodes and optiedges.
Partition(optigraph::OptiGraph,node_membership_vector::Vector{Int64},ref_map::Dict)
Create a partition using optigraph
, node_membership_vector
, and 'refmap'. The `refmap` is a mapping of node_indices to the original optinodes.
Partition(optigraph::OptiGraph,optinode_vectors::Vector{Vector{OptiNode}})
Manually create a partition using optigraph
and a vector of vectors containing sets of optinodes that represent each partition.
Plasmo.apply_partition!
— Functionapply_partition!(optigraph::OptiGraph,partition::Partition)
Create subgraphs in optigraph
using a partition
.
Plasmo.aggregate
— Functionaggregate(graph::OptiGraph)
Aggregate the optigraph graph
into a new optinode. Return an optinode and a dictionary which maps optinode variable and constraint references to the original optigraph.
aggregate(graph::OptiGraph,max_depth::Int64)
Aggregate the optigraph 'graph' into a new aggregated optigraph. Return a newly aggregated optigraph and a dictionary which maps new variables and constraints to the original optigraph. max_depth
determines how many levels of subgraphs remain in the new aggregated optigraph. For example, a max_depth
of 0
signifies there should be no subgraphs in the aggregated optigraph.
Plasmo.aggregate!
— Functionaggregate!(graph::OptiGraph, max_depth::Int64)
Aggregate graph
by converting subgraphs into optinodes. The max_depth
determines how many levels of subgraphs remain in the new aggregated optigraph. For example, a max_depth
of 0
signifies there should be no subgraphs in the aggregated optigraph.
LightGraphs.all_neighbors
— FunctionLightGraphs.all_neighbors(graph::OptiGraph, node::OptiNode)
Retrieve the optinode neighbors of node
in the optigraph graph
. Uses an underlying hypergraph to query for neighbors.
LightGraphs.induced_subgraph
— FunctionLightGraphs.induced_subgraph(graph::OptiGraph, nodes::Vector{OptiNode})
Create an induced subgraph of optigraph given a vector of optinodes.
Plasmo.incident_edges
— Functionincident_edges(hypergraph::HyperGraph,hypernode::HyperNode)
Identify the incident hyperedges to a HyperNode
.
incident_edges(hypergraph::HyperGraph,hypernodes::Vector{HyperNode})
Identify the incident hyperedges to a vector of HyperNode
s.
incident_edges(graph::OptiGraph, nodes::Vector{OptiNode})
Retrieve incident edges to a set of optinodes.
incident_edges(graph::OptiGraph, node::OptiNode)
Retrieve incident edges to a single optinode.
Plasmo.induced_edges
— Functioninduced_edges(hypergraph::HyperGraph,hypernodes::Vector{HyperNode})
Identify the induced hyperedges to a vector of HyperNode
s.
NOTE: This currently does not support hypergraphs with unconnected nodes
induced_edges(graph::OptiGraph, nodes::Vector{OptiNode})
Retrieve induced edges to a set of optinodes.
Plasmo.identify_edges
— Functionidentify_edges(hypergraph::HyperGraph,partitions::Vector{Vector{HyperNode}})
Identify both induced partition edges and cut edges given a partition of HyperNode
vectors.
identify_edges(graph::OptiGraph, node_vectors::Vector{Vector{OptiNode}})
Identify induced edges and edge separators from a vector of optinode partitions.
Plasmo.identify_nodes
— Functionidentify_nodes(hypergraph::HyperGraph,partitions::Vector{Vector{HyperEdge}})
Identify both induced partition nodes and cut nodes given a partition of HyperEdge
vectors.
identify_nodes(graph::OptiGraph, node_vectors::Vector{Vector{OptiEdge}})
Identify induced nodes and node separators from a vector of optiedge partitions.
Plasmo.neighborhood
— Functionneighborhood(g::HyperGraph,nodes::Vector{OptiNode},distance::Int64)
Retrieve the neighborhood within distance
of nodes
. Returns a vector of the original vertices and added vertices
neighborhood(graph::OptiGraph, nodes::Vector{OptiNode}, distance::Int64)::Vector{OptiNode})
Return the optinodes within distance
of the given nodes
in the optigraph graph
.
Plasmo.expand
— Functionexpand(graph::OptiGraph, subgraph::OptiGraph, distance::Int64)
Return a new expanded subgraph given the optigraph graph
and an existing subgraph subgraph
. The returned subgraph contains the expanded neighborhood within distance
of the given subgraph
.
Plasmo.hierarchical_edges
— Functionhierarchical_edges(graph::OptiGraph)::Vector{OptiEdge}
Query the edges in graph
that connect its local nodes to nodes in its subgraphs.
Missing docstring for linking_edges
. Check Documenter's build log for details.
Plotting
PlasmoPlots.layout_plot
— FunctionPlasmoPlots.layout_plot(graph::OptiGraph; node_labels = false, subgraph_colors = false, node_colors = false, linewidth = 2.0,linealpha = 1.0, markersize = 30,labelsize = 20, markercolor = :grey,
layout_options = Dict(:tol => 0.01,:C => 2, :K => 4, :iterations => 2),
plt_options = Dict(:legend => false,:framestyle => :box,:grid => false,:size => (800,800),:axis => nothing),
line_options = Dict(:linecolor => :blue,:linewidth => linewidth,:linealpha => linealpha))
Plot a graph layout of the optigraph graph
. The following keyword arguments can be provided to customize the graph layout.
node_labels = false
: whether to label nodes using the corresponding optinode label attribute.subgraph_colors = false
: whether to color nodes according to their subgraph.node_colors = false
: whether to color nodes. Only active ifsubgraph_colors = false
.linewidth = 2.0
: the linewidth attribute for each edge ingraph
.linealpha = 1.0
: the linealpha attribute for each edge ingraph
.markersize = 30
: the markersize which determines the size of each node ingraph
.labelsize = 20
: the size for each node label. Only active ifnode_labels = true
.markercolor = :grey
: the color for each node.layout_options = Dict(:tol => 0.01,:C => 2, :K => 4, :iterations => 2)
: dictionary with options for the layout algorithm.tol
: permitted distance between a current and calculated co-ordinate.C
,K
: scaling parameters.iterations
: number of iterations used to apply forces.
plt_options = Dict(:legend => false,:framestyle => :box,:grid => false,:size => (800,800),:axis => nothing)
: dictionary with primary plotting options.legend
: whether to include legend, or legend position.framestyle
: style of frame used for plot.size
: size of the resulting plot.axis
: whether to include the axis. The axis typically does not make sense for a graph layout plot.- It is also possible to use various plotting options compatible with
Plots.scatter
from thePlots.jl
package.
line_options = Dict(:linecolor => :blue,:linewidth => linewidth,:linealpha => linealpha)
: line plotting options used to display edges in the graph.linecolor
: color to use for each line.linewidth
: linewidth to use for each edge. Defaults to the above option.linealpha
: linealpha to use for each edge. Default to the above option.
PlasmoPlots.matrix_plot
— FunctionPlasmoPlots.matrix_plot(graph::OptiGraph;node_labels = false,labelsize = 24,subgraph_colors = false,node_colors = false,markersize = 1)
Plot a matrix visualization of the optigraph: graph
. The following keyword arguments can be provided to customize the matrix visual.
node_labels = false
: whether to label nodes using the corresponding optinode label attribute.labelsize
: the size for each node label. Only active ifnode_labels = true
.subgraph_colors = false
: whether to color nodes according to their subgraph.node_colors = false
: whether to color nodes. Only active ifsubgraph_colors = false
.markersize = 1
: Size of the linking constraints in the matrix representation.