STAG C++  1.3.0
Spectral Toolkit of Algorithms for Graphs
Loading...
Searching...
No Matches
stag::Graph Class Reference

#include <graph.h>

Inheritance diagram for stag::Graph:
stag::LocalGraph

Description

The core object used to represent graphs for use with the library.

Graphs are always constructed from sparse matrices, and this is the internal representation used as well. Vertices of the graph are always referred to by their unique integer index. This index corresponds to the position of the vertex in the stored adjacency matrix of the graph.

This class supports graphs with positive edge weights. Self-loops are permitted.

Public Member Functions

 Graph (const SprsMat &matrix)
 
 Graph (std::vector< stag_int > &outerStarts, std::vector< stag_int > &innerIndices, std::vector< double > &values)
 
const SprsMatadjacency () const
 
const SprsMatlaplacian ()
 
const SprsMatnormalised_laplacian ()
 
const SprsMatsignless_laplacian ()
 
const SprsMatnormalised_signless_laplacian ()
 
const SprsMatdegree_matrix ()
 
const SprsMatinverse_degree_matrix ()
 
const SprsMatlazy_random_walk_matrix ()
 
double total_volume ()
 
double average_degree ()
 
stag_int number_of_vertices () const
 
stag_int number_of_edges () const
 
bool has_self_loops () const
 
bool is_connected ()
 
Graph subgraph (std::vector< stag_int > &vertices)
 
Graph disjoint_union (Graph &other)
 
double degree (stag_int v) override
 
stag_int degree_unweighted (stag_int v) override
 
std::vector< edgeneighbors (stag_int v) override
 
std::vector< stag_intneighbors_unweighted (stag_int v) override
 
std::vector< double > degrees (std::vector< stag_int > vertices) override
 
std::vector< stag_intdegrees_unweighted (std::vector< stag_int > vertices) override
 
bool vertex_exists (stag_int v) override
 
- Public Member Functions inherited from stag::LocalGraph
virtual double degree (stag_int v)=0
 
virtual stag_int degree_unweighted (stag_int v)=0
 
virtual std::vector< edgeneighbors (stag_int v)=0
 
virtual std::vector< stag_intneighbors_unweighted (stag_int v)=0
 
virtual std::vector< double > degrees (std::vector< stag_int > vertices)=0
 
virtual std::vector< stag_intdegrees_unweighted (std::vector< stag_int > vertices)=0
 
virtual bool vertex_exists (stag_int v)=0
 
virtual ~LocalGraph ()=default
 

Constructor & Destructor Documentation

◆ Graph() [1/2]

stag::Graph::Graph ( const SprsMat matrix)
explicit

Create a graph from an Eigen matrix.

The provided matrix should correspond either to the adjacency matrix or Laplacian matrix of the graph. STAG will automatically detect whether the provided matrix is an adjacency matrix or a Laplacian matrix.

Example
#include <iostream>
#include <stag/graph.h>
int main() {
// Construct a sparse matrix representing the
// triangle graph adjacency matrix.
stag_int n = 3;
SprsMat adj(n, n);
adj.coeffRef(0, 1) = 1;
adj.coeffRef(0, 2) = 1;
adj.coeffRef(1, 0) = 1;
adj.coeffRef(1, 2) = 1;
adj.coeffRef(2, 0) = 1;
adj.coeffRef(2, 1) = 1;
// Create a new STAG graph
stag::Graph myGraph(adj);
// Display the adjacency matrix of the graph
std::cout << *myGraph.adjacency() << std::endl;
return 0;
}
The core object used to represent graphs for use with the library.
Definition: graph.h:179
long long stag_int
Definition: graph.h:26
Eigen::SparseMatrix< double, Eigen::ColMajor, stag_int > SprsMat
Definition: graph.h:32

The provided matrix must be symmetric, and may include self-loops.

Parameters
matrixthe sparse eigen matrix representing the adjacency matrix or Laplacian matrix of the graph.
Exceptions
domain_errorif the provided matrix is not symmetric

◆ Graph() [2/2]

stag::Graph::Graph ( std::vector< stag_int > &  outerStarts,
std::vector< stag_int > &  innerIndices,
std::vector< double > &  values 
)

Create a graph from raw arrays describing a CSC sparse matrix.

To use this constructor, you should understand the CSC sparse matrix format. Note that this library uses the ColMajor format from the Eigen library. For more information, refer to the Eigen Documentation.

Parameters
outerStartsthe indices of the start of each row in the CSR matrix
innerIndicesthe column indices of each non-zero element in the matrix
valuesthe values of each non-zero element in the matrix

Member Function Documentation

◆ adjacency()

const SprsMat * stag::Graph::adjacency ( ) const

Return the sparse adjacency matrix of the graph.

Returns
a sparse Eigen matrix representing the graph adjacency matrix.

◆ laplacian()

const SprsMat * stag::Graph::laplacian ( )

Return the Laplacian matrix of the graph.

The Laplacian matrix is defined by

\[ L = D - A \]

where \(D\) is the diagonal matrix of vertex degrees (stag::Graph::degree_matrix) and A is the adjacency matrix of the graph (stag::Graph::adjacency).

Returns
a sparse Eigen matrix representing the graph Laplacian

◆ normalised_laplacian()

const SprsMat * stag::Graph::normalised_laplacian ( )

Return the normalised Laplacian matrix of the graph.

The normalised Laplacian matrix is defined by

\[ \mathcal{L} = D^{-1/2} L D^{-1/2} \]

where \(D\) is the diagonal matrix of vertex degrees (stag::Graph::degree_matrix) and \(L\) is the Laplacian matrix of the graph (stag::Graph::laplacian).

Returns
a sparse Eigen matrix representing the normalised Laplacian

◆ signless_laplacian()

const SprsMat * stag::Graph::signless_laplacian ( )

Return the signless Laplacian matrix of the graph.

The signless Laplacian matrix is defined by

\[ J = D + A \]

where \(D\) is the diagonal matrix of vertex degrees (stag::Graph::degree_matrix) and A is the adjacency matrix of the graph (stag::Graph::adjacency).

Returns
a sparse Eigen matrix representing the signless graph Laplacian

◆ normalised_signless_laplacian()

const SprsMat * stag::Graph::normalised_signless_laplacian ( )

Return the normalised signless Laplacian matrix of the graph.

The normalised signless Laplacian matrix is defined by

\[ \mathcal{J} = D^{-1/2} J D^{-1/2} \]

where \(D\) is the diagonal matrix of vertex degrees (stag::Graph::degree_matrix) and \(J\) is the signless Laplacian matrix of the graph (stag::Graph::signless_laplacian).

Returns
a sparse Eigen matrix representing the normalised Laplacian

◆ degree_matrix()

const SprsMat * stag::Graph::degree_matrix ( )

The degree matrix of the graph.

The degree matrix is an \(n \times n\) matrix such that each diagonal entry is the degree of the corresponding node.

Returns
a sparse Eigen matrix

◆ inverse_degree_matrix()

const SprsMat * stag::Graph::inverse_degree_matrix ( )

The inverse degree matrix of the graph.

The inverse degree matrix is an \(n \times n\) matrix such that each diagonal entry is the inverse of the degree of the corresponding node, or 0 if the node has degree 0.

Returns
a sparse Eigen matrix

◆ lazy_random_walk_matrix()

const SprsMat * stag::Graph::lazy_random_walk_matrix ( )

The lazy random walk matrix of the graph.

The lazy random walk matrix is defined to be

\[ \frac{1}{2} I + \frac{1}{2} A D^{-1} \]

where \(I\) is the identity matrix, \(A\) is the graph adjacency matrix and \(D\) is the degree matrix of the graph.

Returns
a sparse Eigen matrix

◆ total_volume()

double stag::Graph::total_volume ( )

The total volume of the graph.

The volume is defined as the sum of the node degrees.

Returns
the graph's volume.

◆ average_degree()

double stag::Graph::average_degree ( )

The average degree of the graph.

This is defined as the sum of the node degrees divided by the number of nodes.

Returns
the graph's average degree.

◆ number_of_vertices()

stag_int stag::Graph::number_of_vertices ( ) const

The number of vertices in the graph.

◆ number_of_edges()

stag_int stag::Graph::number_of_edges ( ) const

The number of edges in the graph.

This is defined based on the number of non-zero elements in the adjacency matrix, and ignores the weights of the edges.

◆ has_self_loops()

bool stag::Graph::has_self_loops ( ) const

Returns a boolean indicating whether this graph contains self loops.

◆ is_connected()

bool stag::Graph::is_connected ( )

Returns a boolean indicating whether the graph is connected.

The running time of this method is \(O(m)\), where \(m\) is the number of edges in the graph.

◆ subgraph()

Graph stag::Graph::subgraph ( std::vector< stag_int > &  vertices)

Construct and return a subgraph of this graph.

Note that the vertex indices will be changed in the subgraph.

Parameters
verticesthe vertices in the induced subgraph
Returns
a new stag::Graph object representing the subgraph induced by the given vertices

◆ disjoint_union()

Graph stag::Graph::disjoint_union ( Graph other)

Construct and return the disjoint union of this graph and another.

The disjoint union of two graphs \(G\) and \(H\) is a graph containing \(G\) and \(H\) as disconnected subgraphs.

Parameters
otherthe other graph to be combined with this one
Returns
a new stag::Graph object representing the union of this graph with the other one

◆ degree()

double stag::Graph::degree ( stag_int  v)
overridevirtual

Given a vertex v, return its weighted degree.

A self-loop of weight \(1\) contributes \(2\) to the vertex's degree.

Implements stag::LocalGraph.

◆ degree_unweighted()

stag_int stag::Graph::degree_unweighted ( stag_int  v)
overridevirtual

Given a vertex v, return its unweighted degree. That is, the number of neighbors of v, ignoring the edge weights.

Implements stag::LocalGraph.

◆ neighbors()

std::vector< edge > stag::Graph::neighbors ( stag_int  v)
overridevirtual

Given a vertex v, return a vector of edges representing the neighborhood of v.

The returned edges will all have the ordering (v, x) such that edge.v = v.

Parameters
van int representing some vertex in the graph
Returns
an edge vector containing the neighborhood information

Implements stag::LocalGraph.

◆ neighbors_unweighted()

std::vector< stag_int > stag::Graph::neighbors_unweighted ( stag_int  v)
overridevirtual

Given a vertex v, return a vector containing the neighbors of v.

The weights of edges to the neighbors are not returned by this method.

Parameters
van int representing some vertex in the graph
Returns
an int vector giving the neighbors of v

Implements stag::LocalGraph.

◆ degrees()

std::vector< double > stag::Graph::degrees ( std::vector< stag_int vertices)
overridevirtual

Given a list of vertices, return the degrees of each vertex in the list.

Providing a method for computing the degrees 'in bulk' increases the efficiency of algorithms on graphs which are not stored directly in memory.

Parameters
verticesa vector of ints representing the vertices to be queried.
Returns
a vector of degrees

Implements stag::LocalGraph.

◆ degrees_unweighted()

std::vector< stag_int > stag::Graph::degrees_unweighted ( std::vector< stag_int vertices)
overridevirtual

Given a list of vertices, return their unweighted degrees.

Parameters
verticesa vector of ints representing the vertices to be queried.
Returns
a vector of integer degrees

Implements stag::LocalGraph.

◆ vertex_exists()

bool stag::Graph::vertex_exists ( stag_int  v)
overridevirtual

Given a vertex ID, returns true or false to indicate whether the vertex exists in the graph.

Parameters
vthe vertex index to check
Returns
a boolean indicating whether there exists a vertex with the given index

Implements stag::LocalGraph.