application: graph

The application graph deals with directed and undirected graphs. They can be defined abstractly as a set of nodes and EDGES or as part of a larger structure for instance as the vertex-edge graph of a polytope.

imports from: common

Objects

  •  
    Category: Combinatorics
    UNDOCUMENTED
    derived from: Graph<Undirected>

    Properties of GeometricGraph

    •  
      BOUNDING_BOX: common::Matrix

      Since a Voronoi polyhedron is unbounded it must be artificially bounded for visualization purposes. Allowed is any set of hyperplanes which makes the projection onto the last d-1 coordinates bounded. By default, these are the vertical facets of a suitably scaled cube.

    •  
      COORDINATES: common::NodeMap

      The coordinates of the nodes of the graph.

  •  
    Category: Combinatorics

    A graph with optional node and edge attributes.

    Specializations of Graph

    Properties of Graph

    •  
      ADJACENCY: common::Graph

      combinatorial description of the Graph in the form of adjacency matrix


      Example:
      • This prints the adjacency matrix of the complete graph of three nodes, whereby the i-th node is connected to the nodes shown in the i-th row:> print complete(3)->ADJACENCY; {1 2} {0 2} {0 1}
    •  
      AVERAGE_DEGREE: common::Rational
      Only defined for Graph<Undirected>

      The average degree of a node


      Example:
      • The following prints the average degree of a node for the complete bipartite graph with 1 and 3 nodes:> print complete_bipartite(1,3)->AVERAGE_DEGREE; 3/2
    •  
      BICONNECTED_COMPONENTS: common::IncidenceMatrix<NonSymmetric>
      Only defined for Graph<Undirected>

      Biconnected components. Every row contains nodes from a single component. Articulation nodes may occur in several rows.


      Example:
      • The following prints the biconnected components of two triangles having a single common node: > $g = new Graph<Undirected>(ADJACENCY=>[[1,2],[0,2],[0,1,3,4],[2,4],[2,3]]);> print $g->BICONNECTED_COMPONENTS; {2 3 4} {0 1 2}
    •  
      BIPARTITE: common::Bool
      Only defined for Graph<Undirected>

      True if the graph is a bipartite.


      Example:
      • The following checks if a square (as a graph) is bipartite:> $g = new Graph<Undirected>(ADJACENCY=>[[1,3],[0,2],[1,3],[0,2]]);> print $g->BIPARTITE; true
    •  
      CHARACTERISTIC_POLYNOMIAL: common::UniPolynomial<Rational, Int>

      Characteristic polynomial of the adjacency matrix; its roots are the eigenvalues


      Example:
      • The following prints the characteristic polynomial of a complete graph with three nodes:> print complete(3)->CHARACTERISTIC_POLYNOMIAL x^3 -3*x -2
    •  
      CONNECTED: common::Bool

      True if the graph is a connected graph.


      Example:
      • The following checks whether a graph is connected or not. A single edge with its endpoints is obviously connected. We construct the adjacency matrix beforehand:> $IM = new IncidenceMatrix<Symmetric>([[1],[0]]);> print new Graph(ADJACENCY=>$IM)->CONNECTED; true The same procedure yields the opposite outcome after adding an isolated node:> $IM = new IncidenceMatrix<Symmetric>([[1],[0],[]]);> print new Graph(ADJACENCY=>$IM)->CONNECTED; false
    •  
      CONNECTED_COMPONENTS: common::IncidenceMatrix<NonSymmetric>
      Only defined for Graph<Undirected>

      The connected components. Every row contains nodes from a single component.


      Example:
      • The following prints the connected components of two seperate triangles:> $g = new Graph<Undirected>(ADJACENCY=>[[1,2],[0,2],[0,1],[4,5],[3,5],[3,4]]);> print $g->CONNECTED_COMPONENTS; {0 1 2} {3 4 5}
    •  
      CONNECTIVITY: common::Int
      Only defined for Graph<Undirected>

      Node connectivity of the graph, that is, the minimal number of nodes to be removed from the graph such that the result is disconnected.


      Example:
      • The following prints the connectivity of the complete bipartite graph of 2 and 4 nodes:> print complete_bipartite(2,4)->CONNECTIVITY; 2
    •  
      DEGREE_SEQUENCE: common::Map<Int, Int>
      Only defined for Graph<Undirected>

      How many times a node of a given degree occurs


      Example:
      • The following prints how often each degree of a node in the complete bipartite graph with 1 and 3 nodes occurs, whereby the first entry of a tuple represents the degree:> print complete_bipartite(1,3)->DEGREE_SEQUENCE; {(1 3) (3 1)}
    •  
      DIAMETER: common::Int

      The diameter of the graph.


      Example:
      • This prints the diameter of the complete bipartite graph of 1 and 3 nodes:> print complete_bipartite(1,3)->NODE_DEGREES; 3 1 1 1
    •  
      EIGENVALUES_LAPLACIAN: common::Vector<Float>

      eigenvalues of the discrete laplacian operator of a graph


      Example:
      • The following prints the eigenvalues of the discrete laplacian operator of the complete graph with 3 nodes:> print complete(3)->EIGENVALUES_LAPLACIAN 3 -3 0
    •  
      MAX_CLIQUES: common::PowerSet<Int>
      Only defined for Graph<Undirected>

      The maximal cliques of the graph, encoded as node sets.


      Example:
      • The following prints the maximal cliques of two complete graphs with 3 nodes being connected by a single edge:> $g = new Graph<Undirected>(ADJACENCY=>[[1,2],[0,2],[0,1,3],[2,4,5],[3,5],[3,4]]);> print $g->MAX_CLIQUES; {{0 1 2} {2 3} {3 4 5}}
    •  
      NODE_DEGREES: common::Array<Int>

      Degrees of nodes in the graph.


      Example:
      • This prints the node degrees of the complete bipartite graph of 1 and 3 nodes:> print complete_bipartite(1,3)->NODE_DEGREES; 3 1 1 1
    •  
      NODE_IN_DEGREES: common::Array<Int>
      Only defined for Graph<Directed>

      The number of incoming edges of the graph nodes.


      Example:
      • To print the number of incoming edges of a directed version of the complete graph with 3 nodes, type this:> $g = new Graph<Directed>(ADJACENCY=>[[1,2],[2],[]]);> print $g->NODE_IN_DEGREES; 0 1 2
    •  
      NODE_LABELS: common::Array<String>

      Labels of the nodes of the graph.


      Example:
      • The following prints the node labels of the generalized_johnson_graph with parameters (4,2,1):> print generalized_johnson_graph(4,2,1)->NODE_LABELS; {0 1} {0 2} {0 3} {1 2} {1 3} {2 3}
    •  
      NODE_OUT_DEGREES: common::Array<Int>
      Only defined for Graph<Directed>

      The number of outgoing edges of the graph nodes.


      Example:
      • To print the number of incoming edges of a directed version of the complete graph with 3 nodes, type this:> $g = new Graph<Directed>(ADJACENCY=>[[1,2],[2],[]]);> print $g->NODE_OUT_DEGREES; 2 1 0
    •  
      N_CONNECTED_COMPONENTS: common::Int
      Only defined for Graph<Undirected>

      Number of CONNECTED_COMPONENTS of the graph.


      Example:
      • The following prints the number of connected components of two seperate triangles:> $g = new Graph<Undirected>(ADJACENCY=>[[1,2],[0,2],[0,1],[4,5],[3,5],[3,4]]);> print $g->N_CONNECTED_COMPONENTS; 2
    •  
      N_EDGES: common::Int

      Number of EDGES of the graph.


      Example:
      • This prints the number of edges of the complete graph of 4 nodes (which is 4 choose 2):> print complete(4)->N_EDGES; 6
    •  
      N_NODES: common::Int

      Number of nodes of the graph.


      Example:
      • This prints the number of nodes of the complete bipartite graph of 1 and 3 nodes (which is the sum of the nodes in each partition):> print complete_bipartite(1,3)->N_NODES; 4
    •  
      SIGNATURE: common::Int
      Only defined for Graph<Undirected>

      Difference of the black and white nodes if the graph is BIPARTITE. Otherwise -1.


      Example:
      • The following prints the signature of the complete bipartite graph with 2 and 7 nodes:> print complete_bipartite(2,7)->SIGNATURE; 5
    •  
      SIGNED_INCIDENCE_MATRIX: common::SparseMatrix<Int, NonSymmetric>

      Signed vertex-edge incidence matrix; for undirected graphs, the orientation comes from the lexicographic order of the nodes.


      Example:
      • The following prints the signed incidence matrix of the cylcic graph with 4 nodes:> print cycle_graph(4)->SIGNED_INCIDENCE_MATRIX; 1 0 1 0 -1 1 0 0 0 -1 0 1 0 0 -1 -1
    •  
      STRONGLY_CONNECTED: common::Bool
      Only defined for Graph<Directed>

      True if the graph is strongly connected


      Example:
      • The following checks a graph for strong connectivity. First we construct a graph with 4 nodes consisting of two directed circles, which are connected by a single directed edge. We then and check the property:> $g = new Graph<Directed>(ADJACENCY=>[[1],[0,2],[3],[2]]);> print $g->STRONGLY_CONNECTED; false The same procedure yields the opposite result for the graph with the reversed edge added in the middle: > $g = new Graph<Directed>(ADJACENCY=>[[1],[0,2],[1,3],[2]]);> print $g->STRONGLY_CONNECTED; true
    •  
      STRONG_COMPONENTS: common::IncidenceMatrix<NonSymmetric>
      Only defined for Graph<Directed>

      Strong components. Every row contains nodes from a single component


      Example:
      • To print the strong connected components of a graph with 4 nodes consisting of two directed circles and which are connected by a single directed edge, type this:> $g = new Graph<Directed>(ADJACENCY=>[[1],[0,2],[3],[2]]);> print $g->STRONG_COMPONENTS; {2 3} {0 1}
    •  
      TRIANGLE_FREE: common::Bool
      Only defined for Graph<Undirected>

      Determine whether the graph has triangles or not.


      Example:
      • The following checks if the petersen graph contains triangles:> print petersen()->TRIANGLE_FREE; true
    •  
      WEAKLY_CONNECTED: common::Bool
      Only defined for Graph<Directed>

      True if the graph is weakly connected


      Example:
      • The following checks a graph for weak connectivity. First we construct a graph with 4 nodes consisting of two directed circles, which are connected by a single directed edge. We then check the property:> $g = new Graph<Directed>(ADJACENCY=>[[1],[0,2],[3],[2]]);> print $g->WEAKLY_CONNECTED; true
    •  
      WEAKLY_CONNECTED_COMPONENTS: common::IncidenceMatrix<NonSymmetric>
      Only defined for Graph<Directed>

      Weakly connected components. Every row contains nodes from a single component


      Example:
      • To print the weakly connected components of a graph with 4 nodes consisting of two directed circles, which are connected by a single directed edge, type this:> $g = new Graph<Directed>(ADJACENCY=>[[1],[0,2],[3],[2]]);> print $g->WEAKLY_CONNECTED_COMPONENTS; {0 1 2 3} The same procedure yields the opposite outcome using the same graph without the linking edge between the two circles:> $g = new Graph<Directed>(ADJACENCY=>[[1],[0],[3],[2]]);> print $g->WEAKLY_CONNECTED_COMPONENTS; {0 1} {2 3}

    User Methods of Graph

    •  

      These methods capture combinatorial information of the object. Combinatorial properties only depend on combinatorial data of the object like, e.g., the face lattice.

      •  
        EDGES () → Array<Set<Int>>

        Explore the graph as a sequence of its edges.

    •  

      These methods are for visualization.

      •  
        VISUAL () → Visual::Graph

        Visualizes the graph. Decorations may include Coord, disabling the default spring embedder.

        Options
        Intseed
        random seed value for the spring embedder
        option list:Visual::Graph::decorations
        Returns
        Visual::Graph

        Example:
        • The following visualizes the petersen graph with default settings:> petersen()->VISUAL; The following shows some modified visualization style of the same graph:> petersen()->VISUAL(NodeColor=>"green",NodeThickness=>6,EdgeColor=>"purple",EdgeThickness=>4);

    Permutations of Graph

  •  
    Category: Combinatorics

    A Lattice is a poset where join and meet exist for any two elements. It is realized as a directed graph, such that an arbitrar decoration is attached to each node. It is assumed that this decoration always contains at least a face (a Set<Int>) and a rank (an Int). Lattices occur in two different flavors: Sequential and Nonsequential. They are sequential, if it is clear that all nodes occur sorted according to rank (forwards or backwards). Otherwise they should be declared nonsequential. Sequential Lattices are more efficient in terms of storing them in XML files.

    derived from: Graph<Directed>

    Specializations of Lattice

    Properties of Lattice

    •  

      These properties capture combinatorial information of the object. Combinatorial properties only depend on combinatorial data of the object like, e.g., the face lattice.

      •  
        BOTTOM_NODE: common::Int

        The index of the bottom node


        Example:
        • The following prints the bottom node of the face lattice of the 2-simplex (triangle):> print simplex(2)->HASSE_DIAGRAM->BOTTOM_NODE; 0
      •  
        DECORATION: common::NodeMap

        This is the data associated to each node. The prototype for this is BasicDecoration, which consists of properties face and rank.


        Example:
        • The following prints this property of the face lattice of the 2-simplex (triangle):> print simplex(2)->HASSE_DIAGRAM->DECORATION; ({} 0) ({0} 1) ({1} 1) ({2} 1) ({1 2} 2) ({0 2} 2) ({0 1} 2) ({0 1 2} 3)
      •  
        DIMS: common::Array<Int>

        Kept only for backwards compatibility. Basically encodes the INVERSE_RANK_MAP in FaceLattice objects prior to 3.0.7

      •  
        FACES: common::NodeMap<Directed, Set<Int>>

        The face of each node, realized as a NodeMap. This property is kept for two reasons: As a convenient way to access only the face part of the decoration (in this case the property is temporary) and for reasons of backwards compatibility.


        Example:
        • The following prints the faces of the face lattice of the 2-simplex (triangle):> print simplex(2)->HASSE_DIAGRAM->FACES; {} {0} {1} {2} {1 2} {0 2} {0 1} {0 1 2}
      •  
        INVERSE_RANK_MAP: InverseRankMap

        This property provides an efficient way to enumerate all nodes of a given rank. Internally these are realized differently, depending on whether the Lattice is Sequential or Nonsequential. Both provide the same user methods though.


        Example:
        • The following prints this property of the face lattice of the 2-simplex (triangle), where the tuples represent the ranges of nodes belonging to a specific rank:> print simplex(2)->HASSE_DIAGRAM->INVERSE_RANK_MAP; {(0 (0 0)) (1 (1 3)) (2 (4 6)) (3 (7 7))}
      •  
        TOP_NODE: common::Int

        The index of the top node


        Example:
        • The following prints the top node of the face lattice of the 2-simplex (triangle):> print simplex(2)->HASSE_DIAGRAM->TOP_NODE; 7

    User Methods of Lattice

    •  

      These methods capture combinatorial information of the object. Combinatorial properties only depend on combinatorial data of the object like, e.g., the face lattice.

      •  
        dual_faces () → Array<Set<Int> >

        UNDOCUMENTED
        Returns
        Array<Set<Int> >
        For each node, contains the indices of maximal nodes it lies below.

        Example:
        • The following prints the dual faces of the face lattice of the 2-simplex (triangle):> print simplex(2)->HASSE_DIAGRAM->dual_faces(); {0 1 2} {1 2} {0 2} {0 1} {0} {1} {2} {}
      •  
        nodes_of_rank (r) → List<Int>

        UNDOCUMENTED
        Parameters
        Intr
        Returns
        List<Int>
        All indices of nodes of rank r

        Example:
        • The following prints the nodes of rank 1 of the face lattice of the 2-simplex (triangle):> print simplex(2)->HASSE_DIAGRAM->nodes_of_rank(1); {1 2 3}
      •  
        nodes_of_rank_range (r1, r2) → List<Int>

        UNDOCUMENTED
        Parameters
        Intr1
        Intr2
        Returns
        List<Int>
        or Set<Int> All indices of rank r1 <= r <= r2

        Example:
        • The following prints the nodes with rank between 1 and 2 of the face lattice of the 2-simplex (triangle):> print simplex(2)->HASSE_DIAGRAM->nodes_of_rank_range(1,2); {1 2 3 4 5 6}
      •  
        rank () → Int

        UNDOCUMENTED
        Returns
        Int
        The rank of the TOP_NODE

        Example:
        • The following prints the rank of the top node of the face lattice of the 2-simplex (triangle):> print simplex(2)->HASSE_DIAGRAM->rank(); 3
    •  

      These methods are for visualization.

      •  
        VISUAL () → Visual::Lattice

        Visualize the Lattice.

        Options
        Intseed
        random seed value for the node placement
        option list:Visual::Lattice::decorations
        Returns
        Visual::Lattice

        Example:
        • The following visualizes the face lattice of the 2-simplex (triangle) with default settings:> simplex(2)->HASSE_DIAGRAM->VISUAL; The following shows some modified visualization style of the same lattice:> simplex(2)->HASSE_DIAGRAM->VISUAL(NodeColor=>"green",EdgeThickness=>2,EdgeColor=>"purple");
      •  
        VISUAL_DUAL () → Visual::Lattice

        Visualize the dual Lattice. This only produces meaningful results for lattice where the codimension one nodes generate the lattice under intersection.

        Options
        Intseed
        random seed value for the node placement
        option list:Visual::Lattice::decorations
        Returns
        Visual::Lattice

        Example:
        • The following visualizes the dual face lattice of the 2-simplex (triangle) with default settings:> simplex(2)->HASSE_DIAGRAM->VISUAL_DUAL; The following shows some modified visualization style of the same lattice:> simplex(2)->HASSE_DIAGRAM->VISUAL_DUAL(NodeColor=>"green",EdgeThickness=>2,EdgeColor=>"purple");
  •  
    Category: Visualization

    Collection of nodes and edges of an abstract graph amended with visual decoration attributes and an optional embedding in 3-d.

    derived from: Visual::Object
  •  
    Category: Visualization

    Collection of nodes (representing faces of a face lattice) and edges (representing the inclusion relation) amended with visual decoration attributes and an optional embedding in 2-d.

    derived from: Visual::Graph

User Functions

  •  
    graph_from_cycles ()

    UNDOCUMENTED

  •  
    shortest_path_dijkstra (G, weights, source, target, if)

    Find the shortest path in a graph

    Parameters
    GraphG
    a graph without parallel edges
    EdgeMapweights
    edge weights
    Intsource
    the source node
    Inttarget
    the target node
    Boolif
    true, perform backward search
  •  

    These functions capture combinatorial information of the object. Combinatorial properties only depend on combinatorial data of the object like, e.g., the face lattice.

    •  
      all_spanningtrees (G) → Array<Set<int>>

      Calculate all spanning trees for a connected graph along the lines of

      Donald E. Knuth:
      The Art of Computer Programming
      Volume 4, Fascicle 4, 24-31, 2006, Pearson Education Inc.
      Parameters
      GraphG
      beeing connected
      Returns
      Array<Set<int>>

      Example:
      • The following prints all spanning trees of the complete graph with 3 nodes, whereby each line represents a single spanning tree as an edge set:> print all_spanningtrees(complete(3)->ADJACENCY); {0 1} {1 2} {0 2}
    •  
      complement_graph (G) → Graph

      Creates the complement graph of a graph.

      Parameters
      GraphG
      Returns
      Graph

      Example:
      • The following prints the adjancency matrix of the complement graph of the star graph with 4 nodes: > $g = new Graph<Undirected>(ADJACENCY=>[[],[0],[0],[0]]);> print complement_graph($g)->ADJACENCY; {} {2 3} {1 3} {1 2}
    •  
      connectivity (graph) → Int

      Compute the CONNECTIVITY of a given graph using the Ford-Fulkerson flow algorithm.

      Parameters
      props::Graph<Undirected>graph
      Returns
      Int

      Example:
      • Compute the connectivity of the vertex-edge graph of the square:> print connectivity(cube(2)->GRAPH->ADJACENCY); 2 This means that at least two nodes or edges need to be removed in order for the resulting graph not to be connected anymore.
    •  
      eigenvalues_laplacian (G) → Vector<Float>

      Compute the eigenvalues of the discrete Laplacian of a graph.

      Parameters
      GraphG
      Returns
      Vector<Float>

      Example:
      • > $v = eigenvalues_laplacian(cycle_graph(4));> print $v; 4 2 2 0
    •  
      eigenvalues_laplacian (G) → Vector<Float>

      Compute the eigenvalues of the discrete Laplacian of a graph.

      Parameters
      GraphG
      Returns
      Vector<Float>

      Example:
      • > $v = eigenvalues_laplacian(cycle_graph(4)->ADJACENCY);> print $v; 4 2 2 0
    •  
      find_lattice_permutation (L1, L2, permutation) → Permutation

      This takes two lattices and checks whether they are isomorphic, possibly after applying a permutation to the faces. This function only compares faces and ranks of nodes to determine isomorphism

      Parameters
      LatticeL1
      A lattice
      LatticeL2
      Another lattice, having the same decoration and sequential type
      Permutationpermutation
      A permutation to be applied to the faces. If empty, the identity permutation is chosen
      Returns
      Permutation
      A permutation on the nodes of the graph, if the lattices are isomorphic. Otherwise an exeption is thrown.
    •  
      graph_homomorphisms (G, H) → Array<Array<Int>>

      Enumerate all homomorphisms (edge-preserving maps) from one graph to another

      Parameters
      GraphG
      GraphH
      Options
      Boolallow_loops
      Should edges of G be allowed to collapse to a loop when mapped to H? Default 0
      Array<Int>prescribed_map
      A vector of length G.nodes() with those images in G that should be fixed. Negative entries will be enumerated over.
      Returns
      Array<Array<Int>>
    •  
      incidence_matrix (G) → SparseMatrix<Int>

      Compute the unsigned vertex-edge incidence matrix of the graph.

      Parameters
      GraphG
      Returns
      SparseMatrix<Int>

      Example:
      • > $I = incidence_matrix(cycle_graph(4));> print $I 1 0 1 0 1 1 0 0 0 1 0 1 0 0 1 1
    •  
      incidence_matrix (G) → SparseMatrix<Int>

      Compute the unsigned vertex-edge incidence matrix of the graph.


      Example:
      • > $I = incidence_matrix(cycle_graph(4)->ADJACENCY);> print $I; 1 0 1 0 1 1 0 0 0 1 0 1 0 0 1 1
    •  
      laplacian (G) → SparseMatrix<Rational>

      Compute the Laplacian matrix of a graph.

      Parameters
      GraphG
      Returns
      SparseMatrix<Rational>

      Example:
      • > $I = laplacian(cycle_graph(4));> print $I; 2 -1 0 -1 -1 2 -1 0 0 -1 2 -1 -1 0 -1 2
    •  
      laplacian (G) → SparseMatrix<Rational>

      Compute the Laplacian matrix of a graph.

      Parameters
      GraphG
      Returns
      SparseMatrix<Rational>

      Example:
      • > $I = laplacian(cycle_graph(4)->ADJACENCY);> print $I; 2 -1 0 -1 -1 2 -1 0 0 -1 2 -1 -1 0 -1 2
    •  
      lattice_of_chains (lattice) → Lattice<BasicDecoration>

      For a given lattice, this computes the lattice of chains from bottom to top node. The result always includes an artificial top node.

      Parameters
      Lattice<Decoration>lattice
      Returns
      Lattice<BasicDecoration>
      Faces are sets of nodes of elements in the original lattice forming a chain, ranks are lenghts of chains

      Example:
      • The following prints all faces with their ranks of the lattice of chains of the face lattice of the 0-simplex (a single point):> print lattice_of_chains(simplex(0)->HASSE_DIAGRAM)->DECORATION; ({-1} 3) ({0 1} 2) ({0} 1) ({1} 1) ({} 0)
    •  
      line_graph (G) → Graph

      Creates the line graph of a graph.

      Parameters
      GraphG
      Returns
      Graph

      Example:
      • The following prints the adjacency matrix of the line graph of the star graph with 4 nodes:> $g = new Graph<Undirected>(ADJACENCY=>[[],[0],[0],[0]]);> print line_graph($g->ADJACENCY); {1 2} {0 2} {0 1}
    •  
      maximal_chains_of_lattice (F) → IncidenceMatrix

      Computes the set of maximal chains of a Lattice object.

      Parameters
      LatticeF
      Options
      Boolignore_bottom_node
      If true, the bottom node is not included in the chains. False by default
      Boolignore_top_node
      If true, the top node is not included in the chains. False by default
      Returns
      IncidenceMatrix
      Each row is a maximal chain, indices refer to nodes of the Lattice

      Example:
      • The following prints all maximal chains of the face lattice of the 1-simplex (an edge):> print maximal_chains_of_lattice(simplex(1)->HASSE_DIAGRAM); {0 1 3} {0 2 3}
    •  
      n_graph_homomorphisms (G, H) → Int

      Count all homomorphisms (edge-preserving maps) from one graph to another. They are in fact enumerated, but only the count is kept track of using constant memory.

      Parameters
      GraphG
      GraphH
      Options
      Boolallow_loops
      Should edges of G be allowed to collapse to a loop when mapped to H? Default 0
      Array<Int>prescribed_map
      A vector of length G.nodes() with those images in G that should be fixed. Negative entries will be enumerated over.
      Returns
      Int
    •  
      signed_incidence_matrix (G) → SparseMatrix<Int>

      Compute the signed vertex-edge incidence matrix of the graph. In case of undirected graphs, the orientation of the edges is induced by the order of the nodes.

      Parameters
      GraphG
      Returns
      SparseMatrix<Int>

      Example:
      • > $I = signed_incidence_matrix(cycle_graph(4));> print $I; 1 0 1 0 -1 1 0 0 0 -1 0 1 0 0 -1 -1
    •  
      signed_incidence_matrix (G) → SparseMatrix<Int>

      Compute the signed vertex-edge incidence matrix of the graph. In case of undirected graphs, the orientation of the edges is induced by the order of the nodes.


      Example:
      • > $I = signed_incidence_matrix(cycle_graph(4)->ADJACENCY);> print $I; 1 0 1 0 -1 1 0 0 0 -1 0 1 0 0 -1 -1
  •  

    These clients are concerned with automorphisms of graphs and determining whether graphs are isomorphic.

    •  
      automorphisms (graph) → Array<Array<Int>>

      Find the automorphism group of the graph.

      Parameters
      props::Graphgraph
      Returns
      Array<Array<Int>>
      each element encodes a node permutation.

      Example:
      • We first create the vertex-edge graph of the square and then print its automorphism group:> $g=new props::Graph(cube(2)->GRAPH->ADJACENCY);> print automorphisms($g); 0 2 1 3 1 0 3 2 These two permutations generate the group of all node permutations that preserve vertex-edge connectivity.
      Depends on: bliss or nauty
    •  
      automorphisms (m) → Array<Pair<Array<Int>,Array<Int>>>

      Find the automorphism group of the non-symmetric incidence matrix.

      Parameters
      IncidenceMatrix<NonSymmetric>m
      Returns
      Array<Pair<Array<Int>,Array<Int>>>
      each element encodes a permutation of its rows (first) and columns (second).

      Example:
      • The group of combinatorial automorphisms of the 3-cube coincides with the group of (bipartite) graph automorphisms of the vertex/facet incidences. To print this group, type this:> print automorphisms(cube(3)->VERTICES_IN_FACETS); (<0 1 4 5 2 3> <0 1 4 5 2 3 6 7>) (<2 3 0 1 4 5> <0 2 1 3 4 6 5 7>) (<1 0 2 3 4 5> <1 0 3 2 5 4 7 6>) This means that the group is generated by three elements, one per line in the output. Each is written as a pair of permutations. The first gives the action on the facets, the second the action on the vertices.
      Depends on: bliss or nauty
    •  
      automorphisms (m) → Array<Array<Int>>

      Find the automorphism group of the symmetric incidence matrix.

      Parameters
      IncidenceMatrix<Symmetric>m
      Returns
      Array<Array<Int>>
      each element encodes a permutation of its rows (=columns).
      Depends on: bliss or nauty
    •  
      canonical_form (g) → props::Graph

      Find a canonical representation of a graph g. Warning: This representation can depend on the extension (bliss/nauty) used, its version and configuration, as well as the hardware!

      Parameters
      props::Graphg
      Returns
      props::Graph
      Depends on: bliss or nauty
    •  
      canonical_hash (g, k) → Int

      Compute a hash for a graph g independent of the node ordering. Warning: This hash can depend on the extension (bliss/nauty) used, its version and configuration, as well as the hardware! Nauty requires an integer key k as input, bliss will ignore the key.

      Parameters
      props::Graphg
      Intk
      a key for the hash computation, default value 2922320
      Returns
      Int
      Depends on: bliss or nauty
    •  
      canonical_hash (M, k) → Int

      Compute a hash for an incidence matrix I independent of the row ordering. Warning: This hash can depend on the extension (bliss/nauty) used, its version and configuration, as well as the hardware! Nauty requires an integer key k as input, bliss will ignore the key.

      Parameters
      IncidenceMatrixM
      Intk
      a key for the hash computation, default value 2922320
      Returns
      Int
      Depends on: bliss or nauty
    •  
      find_node_permutation (graph1, graph2) → Array<Int>

      Find the node permutation mapping graph1 to graph2. If the given graphs are not isomorphic, throw an expection.

      Parameters
      props::Graphgraph1
      props::Graphgraph2
      Returns
      Array<Int>
      Depends on: bliss or nauty
    •  
      find_row_col_permutation (m1, m2) → Pair<Array<Int>,Array<Int>>

      Find the permutations mapping the non-symmetric incidence matrix m1 to m2. If the given matrices are not isomorphic, throw an expection.

      Parameters
      IncidenceMatrix<NonSymmetric>m1
      IncidenceMatrix<NonSymmetric>m2
      Returns
      Pair<Array<Int>,Array<Int>>
      first permutation applies to the rows, second applies to the columns

      Example:
      • > $m1 = new IncidenceMatrix([1,2],[5,3]);> $m2 = new IncidenceMatrix([4,3],[1,5]);> print find_row_col_permutation($m1,$m2); <1 0> <0 1 4 3 5 2>
      Depends on: bliss or nauty
    •  
      isomorphic (IncidenceMatrix1, IncidenceMatrix2) → Bool

      true if IncidenceMatrix1 and IncidenceMatrix2 are isomorphic.

      Parameters
      IncidenceMatrixIncidenceMatrix1
      IncidenceMatrixIncidenceMatrix2
      Returns
      Bool

      Example:
      • Compare the incidence matrices of the 2-dimensional cube and cross polytope:> $I1 = cube(2)->VERTICES_IN_FACETS;> $I2 = cross(2)->VERTICES_IN_FACETS;> print isomorphic($I1,$I2); true
      Depends on: bliss or nauty
    •  
      isomorphic (graph1, graph2) → Bool

      true if graph1 and graph2 are isomorphic.

      Parameters
      props::Graphgraph1
      props::Graphgraph2
      Returns
      Bool

      Example:
      • Compare the vertex-edge graph of the square with the cycle graph on 4 nodes:> $g1 = cube(2)->GRAPH->ADJACENCY;> $g2 = cycle_graph(4)->ADJACENCY;> print isomorphic($g1,$g2); true
      Depends on: bliss or nauty
    •  
      n_automorphisms (graph) → Int

      Find the order of the automorphism group of the graph.

      Parameters
      props::Graphgraph
      Returns
      Int
      the order of the automorphism group

      Example:
      • > print n_automorphisms(cycle_graph(5)->ADJACENCY); 2
      Depends on: bliss or nauty
  •  

    Special purpose functions.

    •  
      edge_lengths (G, coords) → EdgeMap

      Compute the lengths of all edges of a given graph G from the coordinates coords of its nodes.

      Parameters
      props::Graph<Directed>G
      the input graph
      Matrixcoords
      the coordinates of the nodes
      Returns
      EdgeMap

      Example:
      • The following prints the edge length of the complete graph with 3 nodes and edge lengths given by the coordiantes of the standard 2-simplex:> print edge_lengths(complete(3)->ADJACENCY,simplex(2)->VERTICES); 1 1 1.414213562
    •  
      graph_from_edges (edges) → Graph

      Creates a graph from a given list of edges.

      Parameters
      Array<Set<Int>>edges
      Returns
      Graph

      Example:
      • > $g = graph_from_edges([[1,2],[1,3],[1,4]]);> print $g->ADJACENCY; {} {2 3 4} {1} {1} {1}
  •  

    With these clients you can create special examples of graphs, graphs belonging to parameterized families and random graphs.

    •  
      complete (n) → Graph

      Constructs a complete graph on n nodes.

      Parameters
      Intn
      Returns
      Graph

      Example:
      • To print the adjacency representation of the complete graph on 3 nodes, type this:> print complete(3)->ADJACENCY {1 2} {0 2} {0 1}
    •  
      complete_bipartite (k, l) → Graph

      Constructs a complete bipartite graph on k + l nodes.

      Parameters
      Intk
      Intl
      Returns
      Graph

      Example:
      • To print the adjacency representation of a complete bipartite graph with two nodes per partition, type this:> print complete_bipartite(2,2)->ADJACENCY; {2 3} {2 3} {0 1} {0 1}
    •  
      cycle_graph (n) → Graph

      Constructs a cycle graph on n nodes.

      Parameters
      Intn
      Returns
      Graph

      Example:
      • To print the adjacency representation of the cycle graph on four nodes, type this:> $g = cycle_graph(4);> print $g->ADJACENCY; {1 3} {0 2} {1 3} {0 2}
    •  
      generalized_johnson_graph (n, k, i) → Graph

      Create the generalized Johnson graph on parameters (n,k,i). It has one node for each set in \({[n]}\choose{k}\), and an edge between two nodes iff the intersection of the corresponding subsets is of size i.

      Parameters
      Intn
      the size of the ground set
      Intk
      the size of the subsets
      Inti
      the size of the subsets
      Returns
      Graph

      Example:
      • The following prints the adjacency representation of the generalized johnson graph with the parameters 4,2,1:> print generalized_johnson_graph(4,2,1)->ADJACENCY; {1 2 3 4} {0 2 3 5} {0 1 4 5} {0 1 4 5} {0 2 3 5} {1 2 3 4}
    •  
      johnson_graph (n, k) → Graph

      Create the Johnson graph on parameters (n,k). It has one node for each set in \({[n]}\choose{k}\), and an edge between two nodes iff the intersection of the corresponding subsets is of size k-1.

      Parameters
      Intn
      the size of the ground set
      Intk
      the size of the subsets
      Returns
      Graph

      Example:
      • The following prints the adjacency representation of the johnson graph with the parameters 4,3:> print johnson_graph(4,3)->ADJACENCY; {1 2 3} {0 2 3} {0 1 3} {0 1 2}
    •  
      kneser_graph (n, k) → Graph

      Create the Kneser graph on parameters (n,k). It has one node for each set in \({[n]}\choose{k}\), and an edge between two nodes iff the corresponding subsets are disjoint.

      Parameters
      Intn
      the size of the ground set
      Intk
      the size of the subsets
      Returns
      Graph

      Example:
      • The following prints the adjacency representation of the kneser graph with the parameters 3,1:> print kneser_graph(3,1)->ADJACENCY; {1 2} {0 2} {0 1}
    •  
      neighborhood_graph (D, delta) → Graph

      Constructs the neighborhood graph of a point set S given a parameter delta. The set is passed as its so-called "distance matrix", whose (i,j)-entry is the distance between point i and j. This matrix can e.g. be computed using the distance_matrix function. Two vertices will be adjacent if the distance of the corresponding points is less than delta.

      Parameters
      Matrix<Rational>D
      input point cloud distance matrix (can be upper triangular)
      Rationaldelta
      the maximal distance of neighbored vertices
      Returns
      Graph

      Example:
      • The following prints the neighborhood graph of a distance matrix with a limit of 3.3, producing the graph of a square:> $D = new Matrix<Rational>([[0,17/10,21/10,42/10],[0,0,79/10,31/10],[0,0,0,6/10],[0,0,0,0]]);> print neighborhood_graph($D,3.3)->ADJACENCY; {1 2} {0 3} {0 3} {1 2}
    •  
      path_graph (n) → Graph

      Constructs a path graph on n nodes.

      Parameters
      Intn
      Returns
      Graph
    •  
      petersen () → Graph

      Constructs the Petersen graph.

      Returns
      Graph

      Example:
      • The following prints the adjacency matrix of the petersen graph:> print petersen()->N_NODES; 10
    •  
      random_graph (n) → Graph

      Constructs a random graph with n nodes according to the Erdos-Renyi model. Each edge is chosen uniformly with probability p.

      Parameters
      Intn
      Options
      Rationalp
      the probability of an edge occurring; default 1/2
      Booltry_connected
      whether to try to generate a connected graph, default 1
      Intmax_attempts
      If connected is set, specifies how many times to try to make a connected random graph before giving up.
      Intseed
      controls the outcome of the random number generator; fixing a seed number guarantees the same outcome.
      Returns
      Graph

      Example:
      • The following produces a connected graph on 10 nodes using a specific seed for a random graph model, where an edge between two nodes occurs with probabilty 0.1.> $g = random_graph(10,p=>0.1,try_connected=>1,max_attempts=>50,seed=>100000);> print $g->N_EDGES; 9
  •  

    These functions are for visualization.

    •  
      clip_graph (G, V, BB) → GeometricGraph

      Clip a graph with respect to a given bounding box. Used for the visualization of Voronoi diagrams.

      Parameters
      GraphG
      MatrixV
      MatrixBB
      Returns
      GeometricGraph
    •  
      graphviz (vis_obj ...)

      Draw the given graph or face lattice object using graphviz program neato or dot respectively. The output is rendered in PostScript format and fed into a viewer program, if one is configured. If you prefer to produce another output format, please use the File option and call the neato or dot program manually.

      Parameters
      Visual::Objectvis_obj ...
      objects to display
      Options
      StringFile
      "filename" or "AUTO" Store the graph description in a DOT source file without starting the interactive GUI. The .dot suffix is automatically added to the file name.
      Specify AUTO if you want the filename be automatically derived from the drawing title.
      You can also use any expression allowed for the open function, including "-" for terminal output, "&HANDLE" for an already opened file handle, or "| program" for a pipe.

      Example:
      • The following creates a star graph with 4 nodes and visualizes it via graphviz with default options:> $g = new Graph<Undirected>(ADJACENCY=>[[],[0],[0],[0]]);> graphviz($g->VISUAL); The following shows some modified visualization style of the same graph:> $g = new Graph<Undirected>(ADJACENCY=>[[],[0],[0],[0]]);> graphviz($g->VISUAL(NodeColor=>"green",EdgeColor=>"purple",EdgeThickness=>5));
    •  
      hd_embedder (label_width)

      Create an embedding of the Lattice as a layered graph. The embedding algorithm tries to minimize the weighted sum of squares of edge lengths, starting from a random distribution. The weights are relative to the fatness of the layers. The y-space between the layers is constant.

      Parameters
      Arraylabel_width
      estimates (better upper bounds) of the label width of each node. The computed layout guarantees that the distances between the nodes in a layer are at least equal to the widest label in this layer.
      Options
      Booldual
      the node representing the empty face is put on the topmost level
      Floateps
      calculation accuracy.
      Intseed
      effects the initial placement of the nodes.
    •  
      LEDA_graph (G)

      Write a graph in LEDA input format.

      Parameters
      props::GraphG
    •  
      metapost (vis_obj ...)

      Produce a MetaPost input file with given visual objects.

      Parameters
      Visual::Objectvis_obj ...
      objects to display
      Options
      StringFile
      "filename" or "AUTO" The MetaPost description always has to be stored in a file, there is no interactive viewer for this kind of visualization.
      For the file name you can use any expression allowed for the open function, including "-" for terminal output, "&HANDLE" for an already opened file handle, or "| program" for a pipe. Real file names are automatically completed with the .mp suffix if needed.
      The default setting "AUTO" lets the file name be derived from the drawing title. The automatically generated file name is displayed in the verbose mode.

      Example:
      • The following prints a metapost description of the complete graph with 3 nodes in the console:> metapost(complete(3)->VISUAL,File=>"-");
    •  
      spring_embedder (graph)

      Produce a 3-d embedding for the graph using the spring embedding algorithm along the lines of

      Thomas Fruchtermann and Edward Reingold:
      Graph Drawing by Force-directed Placement.
      Software Practice and Experience Vol. 21, 1129-1164 (1992), no. 11.
      Parameters
      props::Graph<Undirected>graph
      to be embedded.
      Options
      affecting the desired picture
      EdgeMapedge_weights
      relative edge lengths. By default the embedding algorithm tries to stretch all edges to the same length.
      Vectorz-ordering
      an objective function provides an additional force along the z-axis, trying to rearrange nodes in the order of the function growth.
      Floatz-factor
      gain coefficient applied to the z-ordering force.
      Intseed
      random seed for initial node placement on a unit sphere.
      calculation fine-tuning
      Floatscale
      enlarges the ideal edge length
      Floatbalance
      changes the balance between the edge contraction and node repulsion forces
      Floatinertion
      affects how the nodes are moved, can be used to restrain oscillations
      Floatviscosity
      idem
      Floateps
      a threshold for point movement between iterations, below that it is considered to stand still
      Intmax-iterations
      hard limit for computational efforts. The algorithm terminates at latest after that many iterations regardless of the convergence achieved so far.

      Example:
      • The following prints a 3-dimensional embedding of the complete graph on 3 nodes using a specific seed and scaled edge lengths:> print spring_embedder(complete(3)->ADJACENCY, scale=>5, seed=>123); 0.9512273649 -10.00210559 10.36309695 10.61947526 1.391783824 -9.666627553 -11.57070263 8.610321763 -0.6964693941

Property Types

  •  

    These property_types capture combinatorial information of the object. Combinatorial properties only depend on combinatorial data of the object like, e.g., the face lattice.

    •  
      BasicDecoration

      This is the prototype of decorations attached to a Lattice. It consists of two properties, face and rank of type Set<Int> and Int, respectively.

    •  
      InverseRankMap

      This provides an efficient way to obtain all nodes of a lattice of given rank.

      User Methods of InverseRankMap

      •  

        These methods capture combinatorial information of the object. Combinatorial properties only depend on combinatorial data of the object like, e.g., the face lattice.

        •  
          get_map () → Map<Int, List<Int> >

          UNDOCUMENTED
          Returns
          Map<Int, List<Int> >
          or Map<Int, Pair<Int,Int> >. An actual map object sorting nodes according to rank. In the nonsequential case, each integer (= rank) is mapped to a list of the corresponding nodes.\ In the sequential case, it is mapped to the first and last index of all nodes of that rank.
        •  
          nodes_of_rank (r) → List<Int>

          UNDOCUMENTED
          Parameters
          Intr
          Returns
          List<Int>
          All nodes of rank r.
        •  
          nodes_of_rank_range (r1, r2) → List<Int>

          UNDOCUMENTED
          Parameters
          Intr1
          Intr2
          Returns
          List<Int>
          or Set<Int> All indices of rank r1 <= r <= r2
    •  
      Nonsequential

      Nonsequential lattices are those where not all nodes are necessarily sorted to rank.

    •  
      Sequential

      Sequential lattices are those where all nodes are sorted to rank.

Common Option Lists