NWGraph API Reference
Concepts
graph
-
template<typename G>
concept graph - #include <nwgraph/graph_concepts.hpp>
Base concept for types fulfilling requirements of a graph.
The
graph
concept requires the following:The type
G
must becopyable
vertex_id_t<G>
must be a valid typenameIf
g
is of typeG
, the expressionnum_vertices(g)
is a customization point object returning a type convertiblestd::ranges::range_difference_t<G>
To enable standard library types to meet the requirements of the
graph
concept, appropriate definitions ofvertex_id_t
andnum_vertices
have been made in thegraph_concepts.hpp
header file.
adjacency_list_graph
-
template<typename G>
concept adjacency_list_graph - #include <nwgraph/graph_concepts.hpp>
Concept for types fulfilling requirements of a graph in adjacency format.
The
adjacency_list_graph
concept requires the following:The type
G
must meet the requirements ofgraph
The type
G
must meet the requirements ofstd::ranges::random_access_range
The inner range of
G
, i.e., the type given bystd::ranges::range_value_t<G>
must meet the requirements ofstd::ranges::forward_range
The type
vertex_id_t<G>
must be convertible tostd::ranges::range_difference_t
of the inner range.If
g
is of type ‘G’ andu
is an element of typevertex_id_t<G>
, theng[u]
is a valid expression returning a type convertible to the inner range ofG
.If
g
is ofgraph
typeG
ande
is of the value type of the inner range, then the customization point objecttarget(g, e)
returns a type convertible tovertex_id_t<G>
To enable standard library types to meet the requirements of the
adjacency_list_graph
concept, appropriate definitions oftarget
have been made in thegraph_concepts.hpp
header file.Example: Standard library components meeting the requirements of
adjacency_list_graph
#include <forward_list> #include <list> #include <tuple> #include <vector> int main() { static_assert(nw::graph::adjacency_list_graph<std::vector<std::forward_list<std::tuple<int, int, double>>>); static_assert(nw::graph::adjacency_list_graph<std::vector<std::list<std::tuple<int, int, double>>>); static_assert(nw::graph::adjacency_list_graph<std::vector<std::vector<std::tuple<int, int, double>>>); } template <typename G> concept adjacency_list_graph = graph<G>; template <adjacency_list_graph Graph, typename score_t = float, typename accum_t = size_t> std::vector<score_t> brandes_bc(const Graph& G, bool normalize = true) {
degree_enumerable_graph
-
template<typename G>
concept degree_enumerable_graph - #include <nwgraph/graph_concepts.hpp>
Concept for types fulfilling requirements of a graph whose neighborhood size can be queried in constant time.
The
degree_enumerable_graph
concept requires the following:The type
G
must meet the requirements oadjacency_list_graph
If
g
is of graph typeG
andu
is of the vertex id type, the customization point objectdegree(g[u])
returns a type convertible to the difference type ofG
. Thedegree
operation is guaranteed to be constant time.
edge_list_graph
-
template<typename G>
concept edge_list_graph - #include <nwgraph/graph_concepts.hpp>
Concept for types fulfilling requirements of a graph in edge list format.
The
edge_list_graph
concept requires the following:The type
G
must meet the requirements ofgraph
The type
G
must meet the requirements ofstd::ranges::forward_range
If
g
is of type ‘G’ ande
is an element of typevertex_id_t<G>
, thensource(g, e)
returns a type convertible to the vertex type of G.If
g
is of type ‘G’ ande
is an element of typevertex_id_t<G>
, thentarget(g, e)
returns a type convertible to the vertex type of G.
To enable standard library types to meet the requirements of the
edge_list_graph
concept, appropriate definitions ofsource
andtarget
have been made in thegraph_concepts.hpp
header file.Example: Standard library components meeting the requirements of
edge_list_graph
#include <forward_list> #include <list> #include <tuple> #include <vector> int main() { static_assert(nw::graph::edge_list_graph<std::forward_list<std::tuple<int, int, double>>); static_assert(nw::graph::edge_list_graph<std::list<std::tuple<int, int, double>>); static_assert(nw::graph::edge_list_graph<std::vector<std::tuple<int, int, double>>); }
Graph Algorithms
Betweenness Centrality
-
template<adjacency_list_graph Graph, typename score_t = float, typename accum_t = size_t>
std::vector<score_t> nw::graph::brandes_bc(const Graph &G, bool normalize = true) Compute exact betweenness centrality using Brandes algorithm [Bra01].
- Template Parameters
Graph – Type of the graph. Must meet requirements of adjacency_list_graph concept.
score_t – Type of the centrality scores computed for each vertex.
accum_t – Type of path counts.
- Parameters
G – Input graph.
normalize – Flag indicating whether to normalize centrality scores relative to largest score.
- Returns
Vector of centrality for each vertex.
-
template<class score_t, class accum_t, adjacency_list_graph Graph, class OuterExecutionPolicy = std::execution::parallel_unsequenced_policy, class InnerExecutionPolicy = std::execution::parallel_unsequenced_policy>
auto nw::graph::brandes_bc(const Graph &graph, const std::vector<typename Graph::vertex_id_type> &sources, int threads, OuterExecutionPolicy &&outer_policy = {}, InnerExecutionPolicy &&inner_policy = {}, bool normalize = true) Parallel approximate betweenness centrality using Brandes algorithm [Bra01]. Rather than using all vertices in the graph to compute paths, the algorithm uses a specified set of root nodes. Parallelization is effected through C++ standard library execution policies.
- Template Parameters
Graph – Type of the graph. Must meet requirements of adjacency_list_graph concept.
score_t – Type of the centrality scores computed for each vertex.
accum_t – Type of path counts.
OuterExecutionPolicy – Parallel execution policy type for outer loop of algorithm. Default: std::execution::parallel_unsequenced_policy
InnerExecutionPolicy – Parallel execution policy type for inner loop of algorithm. Default: std::execution::parallel_unsequenced_policy
- Parameters
G – Input graph.
normalize – Flag indicating whether to normalize centrality scores relative to largest score.
sources – Vector of starting sources.
outer_policy – Outer loop parallel execution policy.
inner_policy – Inner loop parallel execution policy.
threads – Number of threads being used in computation. Used to compute number of bins in computation.
- Returns
Vector of centrality for each vertex.
Breadth-First Search
-
template<adjacency_list_graph Graph>
auto nw::graph::bfs(const Graph &graph, vertex_id_t<Graph> root) Breadth-First Search.
Perform breadth-first search of a graph, starting from a given vertex.
- Template Parameters
Graph – Type of the input graph. Must meet the requirements of the adjacency_list_graph concept.
- Parameters
graph – The graph to be searched.
root – The starting vertex.
- Returns
The parent list.
-
template<adjacency_list_graph OutGraph, adjacency_list_graph InGraph>
auto nw::graph::bfs(const OutGraph &out_graph, const InGraph &in_graph, vertex_id_t<OutGraph> root, int num_bins = 32, int alpha = 15, int beta = 18) Parallel Breadth-First Search.
Perform parallel breadth-first search of a graph, using Beamer’s direction-optimizing algorithm [BAsanovicP12]. The algorithm requires being able to access the out edges of each vertex as well as the in edges of each vertex. These are passed into the graph as two adajacency lists.
- Template Parameters
OutGraph – Type of graph containing out edges. Must meet the requirements of adjacency_list_graph concept.
InGraph – Type of graph containing in edges. Must meet the requirements of adjacency_list_graph concept.
- Parameters
out_graph – The graph to be searched, representing out edges.
in_graph – The transpose of the graph to be searched, representing in edges.
root – The starting vertex.
num_bins – Number of bins.
alpha – Algorithm parameter.
beta – Algorithm parameter.
- Returns
The parent list.
Connected Components
-
template<typename Execution, adjacency_list_graph Graph1, adjacency_list_graph Graph2>
static auto nw::graph::afforest(Execution &exec, const Graph1 &graph, const Graph2 &t_graph, const size_t neighbor_rounds = 2) Afforest algorithm based on Shiloach-Vishkin with subgraph sampling to skip the largest intermediate component. Labels each vertex in the graph with a component identifier.
- Template Parameters
Execution – execution policy type.
Graph1 – Type of input graph. Must meet requirements of adjacency_list_graph concept.
Graph2 – Type of transpose input graph. Must meet requirements of adjacency_list_graph concept.
- Parameters
exec – Parallel execution policy.
graph – Input graph.
t_graph – Transpose of the input graph.
neighbor_rounds – The number of rounds to do neighborhood subgraph sampling.
- Returns
Vector of component labelings.
Jaccard Similarity
Warning
doxygenfunction: Cannot find function “nw::graph::jaccard_similarity_v0” in doxygen xml output for project “NWgraph” from directory: ./_doxygen/xml
Graph Coloring
-
template<adjacency_list_graph Graph>
void nw::graph::jones_plassmann_coloring(Graph &A, std::vector<size_t> &colors) Graph greedy coloring using maximal independent set.
- Template Parameters
Graph – Type of graph. Must meet the requirements of adjacency_list_graph concept.
- Parameters
A – The input graph.
colors – The array of colors of each vertex.
-
template<adjacency_list_graph Graph>
std::tuple<Unordered_map, size_t> nw::graph::k_core(const Graph &A, int k) A sequential k core algorithm to find the k cores in a graph.
- Template Parameters
Graph – Type of the graph. Must meet requirements of adjacency_list_graph concept.
- Parameters
A – Input graph.
k – The value of k in the k core.
- Returns
std::tuple<Unordered_map, size_t> of an Unordered_map (k core), and the number of vertices in k core.
Minimum Spanning Tree
-
template<edge_list_graph EdgeListT, typename Compare>
EdgeListT nw::graph::kruskal(EdgeListT &E, Compare comp) A sequential Kruskal’s algorithm to find a minimum spanning tree of an undirected edge-weighted graph.
- Template Parameters
EdgeListT – the edge_list_graph graph type
Compare – the comparison function type
- Parameters
E – input edge list
comp – comparison function object for sorting the input edge list
- Returns
EdgeListT output edge list of the minimum spanning tree
-
template<edge_list_graph EdgeListT>
EdgeListT nw::graph::kruskal(EdgeListT &E) A wrapper function to avoid pass compare function as an arg.
- Template Parameters
EdgeListT – the edge_list_graph graph type
- Parameters
E – input edge list
- Returns
EdgeListT output edge list of the minimum spanning tree
-
template<adjacency_list_graph Graph, class Distance, class Weight>
std::vector<vertex_id_t<Graph>> nw::graph::prim(const Graph &graph, vertex_id_t<Graph> source, Weight &&weight) Prim’s minimum spanning tree algorithm A greedy algorithm for a weighted undirected graph.
- Template Parameters
Graph – Graph type. Must meet requirements of adjacency_list_graph.
Distance – Type of edge weight.
Weight – A weight function for a given edge, returns a Distance.
- Parameters
graph – Input graph.
source – Starting vertex.
- Returns
std::vector<vertex_id_t<Graph>>, the predecessor list (the MST).
Max Flow
-
template<typename Dict = default_dict, typename flowtype = double, adjacency_list_graph Graph>
flowtype nw::graph::max_flow(const Graph &A, vertex_id_type source, vertex_id_type sink, size_t max_iters = DEFAULT_MAX) A sequetial algorithm to compute the max flow within a graph.
- Template Parameters
Dict – the dictionary type for max flow
flowtype – the type for flow value
Graph – adjacency_list_graph graph type
- Parameters
A – input graph
source – source vertex
sink – sink vertex
max_iters – the maximum number of iterations
- Returns
flowtype the max flow value
Maximal Independent Set
-
template<adjacency_list_graph Graph>
void nw::graph::maximal_independent_set(const Graph &A, std::vector<size_t> &mis) A sequential algorithm to find maximal independent set in a graph.
- Template Parameters
Graph – input adjacency_list_graph type
- Parameters
A – input graph
mis – the vertices in the maximal independent set
-
template<adjacency_list_graph Graph>
void nw::graph::dag_based_mis(Graph &A, std::vector<bool> &mis) Compute maximal independent set, using DAG adaptor.
- Template Parameters
Graph – Type of the input graph. Must meet requirements of adjacency_list_graph concept.
- Parameters
A – The input graph.
mis – (out) Boolean vector indicating whether corresponding vertex is in the maximal independent set.
PageRank
-
template<adjacency_list_graph Graph, typename Real>
void nw::graph::page_rank(const Graph &graph, const std::vector<typename Graph::vertex_id_type> °rees, std::vector<Real> &page_rank, Real damping_factor, Real threshold, size_t max_iters, size_t num_threads) Parallel page rank.
- Template Parameters
Graph – adjacency_list_graph graph type
Real – page rank score type
- Parameters
graph – input graph
degrees – degree distribution of all vertices
page_rank – container for page rank scores
damping_factor – the probability that an imaginary surfer stops clicking
threshold – error threshold to control converge rate
max_iters – maximum number of iterations to converge
num_threads – number of threads
Single-Source Shortest Paths
-
template<typename Distance, adjacency_list_graph Graph>
std::vector<Distance> nw::graph::dijkstra_er(const Graph &graph, vertex_id_t<Graph> source) Basic Dijkstra’s single-source shortest-paths algorithm, based on bfs_edge_range adaptor.
- Template Parameters
Type – of the edge weights (distances).
Graph – Type of the input graph. Must meet the requirements of the adjacency_list_graph concept.
- Parameters
graph – The input graph.
source – The starting vertex.
- Returns
Vector of distances from the starting node for each vertex in the graph.
-
template<typename Distance, adjacency_list_graph Graph, std::invocable<inner_value_t<Graph>> Weight = std::function<std::tuple_element_t<1, inner_value_t<Graph>>(const inner_value_t<Graph>&)>>
auto nw::graph::dijkstra(const Graph &graph, vertex_id_t<Graph> source, Weight weight =[](auto &e) { return std::get< 1 >(e);}
) Dijkstra’s single-source shortest-paths algorithm, lifted per Lifting Edge Weight.
- Template Parameters
Type – of the edge weights (distances).
Graph – Type of the input graph. Must meet the requirements of the adjacency_list_graph concept.
Weight – Type of function used to compute edge weights.
- Parameters
graph – The input graph.
source – The starting vertex.
weight – Function for computing edge weight.
- Returns
Vector of distances from the starting node for each vertex in the graph.
-
template<class distance_t, adjacency_list_graph Graph, class T, class Weight>
auto nw::graph::delta_stepping(const Graph &graph, vertex_id_t<Graph> source, T delta, Weight weight =[](auto &e) -> auto &{ return std::get< 1 >(e);}
) Delta-stepping single-source shortest-paths.
Sequential implementation of delta-stepping single-source shortest-paths
[MS03].
- Template Parameters
distance_t – Type of distance measure.
Graph – Type of input graph. Must meet the requirements of adjacency_list_graph.
T – Type of delta parameter.
Weight – Type of function used to compute edge weights.
- Parameters
graph – The input graph.
source – The starting vertex.
delta – The delta parameter for the algorithm.
weight – Function to compute weight of an edge.
-
template<class distance_t, adjacency_list_graph Graph, class T>
auto nw::graph::delta_stepping(const Graph &graph, vertex_id_t<Graph> source, T delta) Delta-stepping single-source shortest-paths.
Parallel implementation of delta-stepping single-source shortest-paths [MS03]. Uses Intel TBB for parallelization.
- Template Parameters
distance_t – Type of distance measure.
Graph – Type of input graph. Must meet the requirements of adjacency_list_graph.
T – Type of delta parameter.
- Parameters
graph – The input graph.
source – The starting vertex.
delta – The delta parameter for the algorithm.
Sparse Matrix Sparse Matrix Product
-
template<typename ScalarT, adjacency_list_graph LGraphT, adjacency_list_graph RGraphT, typename MapOpT = std::multiplies<ScalarT>, typename ReduceOpT = std::plus<ScalarT>>
edge_list<directedness::directed, ScalarT> nw::graph::spMatspMat(const LGraphT &A, const RGraphT &B) SpGEMM for C = A * B.
- Todo:
cannot currently pass “const &” for A or B
Need to discuss interface options
- Template Parameters
ScalarT – scalar type
LGraphT – adjacency_list_graph type
RGraphT – adjacency_list_graph type
MapOpT – map operation type
ReduceOpT – reduce operation type
- Parameters
A – Input matrix A
B – Input matrix B
- Returns
edge_list<directedness::directed, ScalarT> a weighted edge list
Triangle Counting
-
template<adjacency_list_graph GraphT>
size_t nw::graph::triangle_count(const GraphT &A) Sequential 2D triangle counting algorithm set intersection the neighbor lists of each pair of vertices.
- Template Parameters
GraphT – adjacency_list_graph
- Parameters
A – graph
- Returns
size_t the number of triangles
-
template<class Op>
std::size_t nw::graph::triangle_count_async(std::size_t threads, Op &&op) Parallel triangle counting using
std::async
.This version of triangle counting uses
threads
std::async
launches to evaluate the passedop
in parallel. Theop
will be provided the thread id, but should capture any other information required to perform the decomposed work.- Template Parameters
Op – The type of the decomposed work.
- Parameters
op – The decomposed work for each
std::async
.- Returns
The += reduced total of counted triangles.
-
template<adjacency_list_graph Graph>
std::size_t nw::graph::triangle_count(const Graph &G, std::size_t threads) Two-dimensional triangle counting, parallel version.
- Template Parameters
Graph – adjacency_list_graph
- Parameters
G – graph
threads – number of threads
- Returns
std::size_t number of triangles
Graph Data Structures
Adjacency List
-
template<int idx, std::unsigned_integral index_type, std::unsigned_integral vertex_id, typename ...Attributes>
class index_adjacency : public nw::graph::unipartite_graph_base, public nw::graph::indexed_struct_of_arrays<index_type, vertex_id, Attributes...> Index adjacency structure. This data structures stores unipartite graph in Compressed Sparse Row format. The underlying data structure is a structure of arrays for storage. Index_adjacency is an adjacency list represenation of a unipartite graph. A unipartite graph has one vertex set.
- Template Parameters
idx – The index type to indicate the type of the graph, can be either 0 or 1.
index_type – The data type used to represent a vertex index, required to be os unsigned integral type.
vertex_id – The data type used to represent a vertex ID, required to be os unsigned integral type.
Attributes – A variadic list of edge property types.
-
template<int idx, typename ...Attributes>
using nw::graph::adjacency = index_adjacency<idx, default_index_t, default_vertex_id_type, Attributes...>
-
template<int idx, typename ...Attributes>
using nw::graph::biadjacency = index_biadjacency<idx, default_index_t, default_vertex_id_type, Attributes...>
Edge List
-
template<std::unsigned_integral vertex_id, typename graph_base_t, directedness direct = directedness::undirected, typename ...Attributes>
class index_edge_list : public graph_base_t, public nw::graph::struct_of_arrays<vertex_id, vertex_id, Attributes...> Index edge list structure. This variadic data structure stores edges with their properties, using a structure of arrays for storage.
An edge in this case is notionally a tuple (u, v, p0, p1, …), where u and v are vertex indices and p0, p1, … are edge properties. If the * directedness tag is set to
directed
, the edges are directed. If the directedness * tag is set toundirected
, the edges are undirected. In either case, a given edge * (u, v, …) is only stored once in the edge list.- Template Parameters
vertex_id – The data type used to represent a vertex index, required to be os unsigned integral type.
graph_base_t – The class of graph represented by the edge list. May be unipartite or bipartite.
directed – The directness of the tag. May be directed or undirected.
Attributes – A variadic list of edge property types.
-
template<directedness edge_directedness = directedness::undirected, typename ...Attributes>
using nw::graph::edge_list = index_edge_list<default_vertex_id_type, unipartite_graph_base, edge_directedness, Attributes...> Type alias for unipartite index edge list structure, using the default_vertex_id_type as the vertex_id type. See
index_edge_list
.
-
template<directedness edge_directedness = directedness::directed, typename ...Attributes>
using nw::graph::bi_edge_list = index_edge_list<default_vertex_id_type, bipartite_graph_base, edge_directedness, Attributes...> Type alias for bipartite index edge list structure, using the default_vertex_id_type as the vertex_id type.
Graph Construction
-
template<edge_list_graph edge_list_t, adjacency_list_graph adj_list_t>
auto nw::graph::fill_adj_list(edge_list_t &el, adj_list_t &al)
-
template<int idx, edge_list_graph edge_list_t, adjacency_list_graph adjacency_t, class Int, class ExecutionPolicy = default_execution_policy>
void nw::graph::fill_directed(edge_list_t &el, Int N, adjacency_t &cs, ExecutionPolicy &&policy = {}) This function fills an adjacency list graph structure from edges contained in a directed edge list.
If edge (u,v) is the edge list, the edge (u,v) will be inserted into the adjacency list — that is, v will be inserted into the neighborhood of u.
If the edges in the edge list have properties, those properties will be copied to the adjacency list.
The adjacency list graph can either be a unipartite or bipartite graph.
- Template Parameters
idx – Which end point to fill in the edge list.
edge_list_t – The type of the edge list.
adjacency_t – The type of the adjacency list.
Int – The type of number of vertices.
ExecutionPolicy – The type of the execution policy.
- Parameters
el – The edge list.
N – Number of vertices at [idx] partition.
cs – The adjacency list.
policy – The parallel execution policy to use in calls to standard library algorithms called by this function.
-
template<int idx, edge_list_graph edge_list_t, class Int, adjacency_list_graph adjacency_t, class ExecutionPolicy = default_execution_policy>
void nw::graph::fill_undirected(edge_list_t &el, Int N, adjacency_t &cs, ExecutionPolicy &&policy = {})
-
template<int idx, edge_list_graph edge_list_t, adjacency_list_graph adjacency_t, class ExecutionPolicy = default_execution_policy>
auto nw::graph::fill(edge_list_t &el, adjacency_t &cs, bool sort_adjacency = false, ExecutionPolicy &&policy = {})
-
template<int idx, edge_list_graph edge_list_t, adjacency_list_graph adjacency_t, class ExecutionPolicy = default_execution_policy>
auto nw::graph::fill(edge_list_t &el, adjacency_t &cs, directedness dir, bool sort_adjacency = false, ExecutionPolicy &&policy = {})
Warning
doxygenfunction: Cannot find function “nw::graph::relabel_by_degree< edge_list_graph edge_list_t, class Vector >” in doxygen xml output for project “NWgraph” from directory: ./_doxygen/xml
Warning
doxygenfunction: Unable to resolve function “nw::graph::relabel_by_degree” with arguments None in doxygen xml output for project “NWgraph” from directory: ./_doxygen/xml. Potential matches:
- template<edge_list_graph edge_list_t, class Vector = std::vector<int>> void relabel_by_degree(edge_list_t &el, std::string direction = "ascending", const Vector °ree = std::vector<int>(0))
- template<int idx, edge_list_graph edge_list_t, class Vector = std::vector<int>> auto relabel_by_degree(edge_list_t &el, std::string direction = "ascending", const Vector °ree = std::vector<int>(0))
-
template<std::ranges::random_access_range V, std::ranges::random_access_range E, adjacency_list_graph Graph = std::vector<std::vector<size_t>>>
auto nw::graph::make_plain_graph(const V &vertices, const E &edges, bool directed = true, size_t idx = 0) Make a plain graph from data, e.g., vector<vector<index>>
-
template<std::ranges::random_access_range V, std::ranges::random_access_range E, adjacency_list_graph Graph = std::vector<std::vector<std::tuple<size_t, size_t>>>>
auto nw::graph::make_index_graph(const V &vertices, const E &edges, bool directed = true, size_t idx = 0) Make an index graph from data, e.g., vector<vector<tuple<index, index>>>
-
template<std::ranges::random_access_range V, std::ranges::forward_range E, adjacency_list_graph Graph = std::vector<std::vector<decltype(std::tuple_cat(std::make_tuple(size_t{}), props(*(begin(E{})))))>>>
auto nw::graph::make_property_graph(const V &vertices, const E &edges, bool directed = true, size_t idx = 0) Make a property graph from data, e.g., vector<vector<tuple<index, properties…>>>
Range Adaptors
-
template<typename Graph>
class back_edge_range
-
template<typename Graph, typename Queue = std::queue<vertex_id_t<Graph>>>
class topdown_bfs_range
-
template<typename Graph>
class bottomup_bfs_range
-
template<typename Graph, typename Queue = std::queue<vertex_id_t<Graph>>>
class bfs_edge_range
-
template<class Iterator>
class cyclic_neighbor_range The cyclic neighbor range adapter.
This adapter takes an underlying random access range and provides the ability to split the range into cycles for TBB. A cycle is a subset of the range such that each subsequent element is some stride from the previous element.
The cyclic range adapter is implemented recursively, that is that each time the range is split it simply returns two ranges that cover the previous range, each with twice the cycle and one offset by one element.
Key to the adapter is the cutoff, which is defined in terms of the maximum stride rather than in terms of the number of elements. A cutoff of
1
implies that the range can’t be split, while a_cutoff_
ofn
means that the range can be split into up ton
cycles.- Template Parameters
Iterator – The type of the underlying range iterator (must be at least random access).
-
template<class Iterator>
class cyclic_range_adaptor stride rather than in terms of the number of elements. A cutoff of
1
implies that the range can’t be split, while a_cutoff_
ofn
means that the range can be split into up ton
cycles.- Template Parameters
Iterator – The type of the underlying range iterator (must be at least random access).
-
template<typename Graph, typename Queue = std::queue<vertex_id_t<Graph>>>
class dag_range
-
template<class Graph, std::size_t... Is>
class edge_range
-
template<typename Graph, typename Queue = std::queue<vertex_id_t<Graph>>, typename Filter = std::function<bool()>>
class filtered_bfs_edge_range
-
template<class Graph, std::size_t... Is>
class neighbor_range
-
template<class Graph, std::size_t... Is>
class plain_range
-
template<typename path_edge, std::integral vertex_id>
class reverse_path
-
template<class Iterator>
class splittable_range_adaptor : public std::ranges::view_base Splittable_range_adaptor can split a range into sub_views of a range.
- Template Parameters
Iterator – iterator type of the range passed in
-
template<class Graph>
class vertex_range Vertex_range assumes the vertex id space is from 0 to n - 1. n is the number of vertices.
Graph Generators
Graph I/O
Warning
doxygenfunction: Unable to resolve function “read_mm” with arguments (const std::string&) in doxygen xml output for project “NWgraph” from directory: ./_doxygen/xml. Potential matches:
- template<directedness sym, typename ...Attributes, edge_list_graph edge_list_t> edge_list_t read_mm(const std::string &filename)
- template<directedness sym, typename ...Attributes> edge_list<sym, Attributes...> read_mm(const std::string &filename)
- template<directedness sym, typename ...Attributes> edge_list<sym, Attributes...> read_mm(std::istream &inputStream)
Warning
doxygenfunction: Unable to resolve function “read_mm” with arguments (const std::string_&) in doxygen xml output for project “NWgraph” from directory: ./_doxygen/xml. Potential matches:
- template<directedness sym, typename ...Attributes, edge_list_graph edge_list_t> edge_list_t read_mm(const std::string &filename)
- template<directedness sym, typename ...Attributes> edge_list<sym, Attributes...> read_mm(const std::string &filename)
- template<directedness sym, typename ...Attributes> edge_list<sym, Attributes...> read_mm(std::istream &inputStream)
-
template<directedness sym, typename ...Attributes>
edge_list<sym, Attributes...> nw::graph::read_mm(std::istream &inputStream)
-
template<size_t w_idx = 0, typename idxtype = void, directedness sym, typename ...Attributes>
void nw::graph::write_mm(const std::string &filename, edge_list<sym, Attributes...> &A, const std::string &file_symmetry = "general")
-
template<size_t w_idx = 0, typename idxtype = void, int idx, typename ...Attributes>
void nw::graph::write_mm(const std::string &filename, adjacency<idx, Attributes...> &A, const std::string &file_symmetry = "general")
-
template<size_t w_idx = 0, typename idxtype = void, int idx, typename ...Attributes>
void nw::graph::write_mm(const std::string &filename, biadjacency<idx, Attributes...> &A, const std::string &file_symmetry = "general")
Containers
-
template<typename ...Attributes>
class array_of_structs : public std::vector<std::tuple<Attributes...>>
-
template<typename index_t, typename ...Attributes>
class indexed_struct_of_arrays Subclassed by nw::graph::index_adjacency< idx, index_type, vertex_id, Attributes >, nw::graph::index_biadjacency< idx, index_type, vertex_id, Attributes >
-
template<std::ranges::random_access_range... Ranges>
struct zipped : public std::tuple<Ranges&...>
Utilities
-
template<typename T = std::size_t>
class counting_output_iterator
-
class par_counting_output_iterator : public std::iterator<std::random_access_iterator_tag, size_t>
-
template<class D = std::chrono::microseconds>
class timer Subclassed by nw::util::life_timer
-
template<typename ThingToSort, typename Comparator, typename IntT, class ExecutionPolicy = std::execution::parallel_unsequenced_policy>
void nw::util::proxysort(const ThingToSort &x, std::vector<IntT> &perm, Comparator comp = std::less<IntT>(), ExecutionPolicy policy = {})
-
template<typename IntT = uint32_t, typename Comparator, typename ThingToSort, class ExecutionPolicy = std::execution::parallel_unsequenced_policy>
auto nw::util::proxysort(const ThingToSort &x, Comparator comp = std::less<IntT>(), ExecutionPolicy policy = {})
Warning
doxygenfunction: Unable to resolve function “nw::graph::intersection_size” with arguments None in doxygen xml output for project “NWgraph” from directory: ./_doxygen/xml. Potential matches:
- template<class A, class B, class C, class D, class ExecutionPolicy> std::size_t intersection_size(A i, B &&ie, C j, D &&je, ExecutionPolicy &&ep)
- template<class A, class B, class C, class D, std::enable_if_t<!std::is_execution_policy_v<std::decay_t<D>>, void**> = nullptr> std::size_t intersection_size(A &&i, B &&ie, C &&j, D &&je)
- template<class A, class B, class Range, class ExecutionPolicy, std::enable_if_t<std::is_execution_policy_v<std::decay_t<ExecutionPolicy>>, void**> = nullptr> std::size_t intersection_size(A &&i, B &&ie, Range &&j, ExecutionPolicy &&ep)
- template<class A, class B, class Range, std::enable_if_t<!std::is_execution_policy_v<std::decay_t<Range>>, void**> = nullptr> std::size_t intersection_size(A &&i, B &&ie, Range &&j)
- template<class R, class S, class ExecutionPolicy, std::enable_if_t<std::is_execution_policy_v<std::decay_t<ExecutionPolicy>>, void**> = nullptr> std::size_t intersection_size(R &&i, S &&j, ExecutionPolicy &&ep)
- template<class R, class S> std::size_t intersection_size(R &&i, S &&j)