dolfin/graph

Documentation for C++ code found in dolfin/graph/*.h

Type definitions

Graph

C++ documentation for Graph from dolfin/graph/Graph.h:

type dolfin::Graph

Vector of unordered Sets.

graph_set_type

C++ documentation for graph_set_type from dolfin/graph/Graph.h:

type dolfin::graph_set_type

Typedefs for simple graph data structures. DOLFIN container for graphs

Classes

BoostGraphColoring

C++ documentation for BoostGraphColoring from dolfin/graph/BoostGraphColoring.h:

class dolfin::BoostGraphColoring

This class colors a graph using the Boost Graph Library.

std::size_t dolfin::BoostGraphColoring::compute_local_vertex_coloring(const Graph &graph, std::vector<ColorType> &colors)

Compute vertex colors.

Parameters:
  • graph
  • colors
std::size_t dolfin::BoostGraphColoring::compute_local_vertex_coloring(const T &graph, std::vector<ColorType> &colors)

Compute vertex colors.

Parameters:
  • graph
  • colors

BoostGraphOrdering

C++ documentation for BoostGraphOrdering from dolfin/graph/BoostGraphOrdering.h:

class dolfin::BoostGraphOrdering

This class computes graph re-orderings. It uses Boost Graph.

T dolfin::BoostGraphOrdering::build_csr_directed_graph(const X &graph)
Parameters:graph
T dolfin::BoostGraphOrdering::build_directed_graph(const X &graph)
Parameters:graph
T dolfin::BoostGraphOrdering::build_undirected_graph(const X &graph)
Parameters:graph
std::vector<int> dolfin::BoostGraphOrdering::compute_cuthill_mckee(const Graph &graph, bool reverse = false)

Compute re-ordering (map[old] -> new) using Cuthill-McKee algorithm

Parameters:
  • graph
  • reverse
std::vector<int> dolfin::BoostGraphOrdering::compute_cuthill_mckee(const std::set<std::pair<std::size_t, std::size_t>> &edges, std::size_t size, bool reverse = false)

Compute re-ordering (map[old] -> new) using Cuthill-McKee algorithm

Parameters:
  • edges
  • size
  • reverse

CSRGraph

C++ documentation for CSRGraph from dolfin/graph/CSRGraph.h:

class dolfin::CSRGraph

Compressed Sparse Row graph. This class provides a Compressed Sparse Row Graph defined by a vector containing edges for each node and a vector of offsets into the edge vector for each node In parallel, all nodes must be numbered from zero on process zero continuously through increasing rank processes. Edges must be defined in terms of the global node numbers. The global node offset of each process is given by node_distribution The format of the nodes, edges and distribution is identical with the formats for ParMETIS and PT-SCOTCH. See the manuals for these libraries for further information.

dolfin::CSRGraph::CSRGraph() = delete

Empty CSR Graph.

dolfin::CSRGraph::CSRGraph(MPI_Comm mpi_comm, const T *xadj, const T *adjncy, std::size_t n)

Create a CSR Graph from ParMETIS style adjacency lists.

Parameters:
  • mpi_comm
  • xadj
  • adjncy
  • n
dolfin::CSRGraph::CSRGraph(MPI_Comm mpi_comm, const std::vector<X> &graph)

Create a CSR Graph from a collection of edges (X is a container some type, e.g. std::vector<unsigned int> or std::set<std::size_t>

Parameters:
  • mpi_comm
  • graph
void dolfin::CSRGraph::calculate_node_distribution()
std::vector<T> &dolfin::CSRGraph::edges()

Vector containing all edges for all local nodes (non-const) (“adjncy” in ParMETIS )

const std::vector<T> &dolfin::CSRGraph::edges() const

Vector containing all edges for all local nodes (“adjncy” in ParMETIS )

std::vector<T> &dolfin::CSRGraph::node_distribution()

Return number of nodes (offset) on each process (non-const)

const std::vector<T> &dolfin::CSRGraph::node_distribution() const

Return number of nodes (offset) on each process.

std::vector<T> &dolfin::CSRGraph::nodes()

Vector containing index offsets into edges for all local nodes (plus extra entry marking end) (“xadj” in ParMETIS )

const std::vector<T> &dolfin::CSRGraph::nodes() const

Vector containing index offsets into edges for all local nodes (plus extra entry marking end) (“xadj” in ParMETIS )

std::size_t dolfin::CSRGraph::num_edges() const

Number of local edges in graph.

std::size_t dolfin::CSRGraph::num_edges(std::size_t i) const

Number of edges from node i.

Parameters:i
const node dolfin::CSRGraph::operator[](std::size_t i) const

Return CSRGraph::node object which provides begin() and end() iterators, also size() , and random-access for the edges of node i.

Parameters:i
std::size_t dolfin::CSRGraph::size() const

Number of local nodes in graph.

T dolfin::CSRGraph::size_global() const

Total (global) number of nodes in parallel graph.

dolfin::CSRGraph::~CSRGraph()

Destructor.

node

C++ documentation for node from dolfin/graph/CSRGraph.h:

class dolfin::CSRGraph::node

Access edges individually by using operator()[] to get a node object

std::vector<T>::const_iterator dolfin::CSRGraph::node::begin() const

Iterator pointing to beginning of edges.

std::vector<T>::const_iterator dolfin::CSRGraph::node::begin_edge
std::vector<T>::const_iterator dolfin::CSRGraph::node::end() const

Iterator pointing to beyond end of edges.

std::vector<T>::const_iterator dolfin::CSRGraph::node::end_edge
dolfin::CSRGraph::node::node(const typename std::vector<T>::const_iterator &begin_it, const typename std::vector<T>::const_iterator &end_it)

Node object, listing a set of outgoing edges.

Parameters:
  • begin_it
  • end_it
const T &dolfin::CSRGraph::node::operator[](std::size_t i) const

Access outgoing edge i of this node.

Parameters:i
std::size_t dolfin::CSRGraph::node::size() const

Number of outgoing edges for this node.

GraphBuilder

C++ documentation for GraphBuilder from dolfin/graph/GraphBuilder.h:

class dolfin::GraphBuilder

This class builds a Graph corresponding to various objects.

Friends: MeshPartitioning.

type dolfin::GraphBuilder::FacetCellMap
std::pair<std::int32_t, std::int32_t> dolfin::GraphBuilder::compute_dual_graph(const MPI_Comm mpi_comm, const boost::multi_array<std::int64_t, 2> &cell_vertices, const CellType &cell_type, const std::int64_t num_global_vertices, std::vector<std::vector<std::size_t>> &local_graph, std::set<std::int64_t> &ghost_vertices)

Build distributed dual graph (cell-cell connections) from minimal mesh data, and return (num local edges, num non-local edges)

Parameters:
  • mpi_comm
  • cell_vertices
  • cell_type
  • num_global_vertices
  • local_graph
  • ghost_vertices
std::int32_t dolfin::GraphBuilder::compute_local_dual_graph(const MPI_Comm mpi_comm, const boost::multi_array<std::int64_t, 2> &cell_vertices, const CellType &cell_type, std::vector<std::vector<std::size_t>> &local_graph, FacetCellMap &facet_cell_map)
Parameters:
  • mpi_comm
  • cell_vertices
  • cell_type
  • local_graph
  • facet_cell_map
std::int32_t dolfin::GraphBuilder::compute_local_dual_graph_keyed(const MPI_Comm mpi_comm, const boost::multi_array<std::int64_t, 2> &cell_vertices, const CellType &cell_type, std::vector<std::vector<std::size_t>> &local_graph, FacetCellMap &facet_cell_map)
Parameters:
  • mpi_comm
  • cell_vertices
  • cell_type
  • local_graph
  • facet_cell_map
std::int32_t dolfin::GraphBuilder::compute_nonlocal_dual_graph(const MPI_Comm mpi_comm, const boost::multi_array<std::int64_t, 2> &cell_vertices, const CellType &cell_type, const std::int64_t num_global_vertices, std::vector<std::vector<std::size_t>> &local_graph, FacetCellMap &facet_cell_map, std::set<std::int64_t> &ghost_vertices)
Parameters:
  • mpi_comm
  • cell_vertices
  • cell_type
  • num_global_vertices
  • local_graph
  • facet_cell_map
  • ghost_vertices
Graph dolfin::GraphBuilder::local_graph(const Mesh &mesh, const GenericDofMap &dofmap0, const GenericDofMap &dofmap1)

Build local graph from dofmap.

Parameters:
  • mesh
  • dofmap0
  • dofmap1
Graph dolfin::GraphBuilder::local_graph(const Mesh &mesh, const std::vector<std::size_t> &coloring_type)

Build local graph from mesh (general version)

Parameters:
  • mesh
  • coloring_type
Graph dolfin::GraphBuilder::local_graph(const Mesh &mesh, std::size_t dim0, std::size_t dim1)

Build local graph (specialized version)

Parameters:
  • mesh
  • dim0
  • dim1

GraphColoring

C++ documentation for GraphColoring from dolfin/graph/GraphColoring.h:

class dolfin::GraphColoring

This class provides a common interface to graph coloring libraries.

std::size_t dolfin::GraphColoring::compute_local_vertex_coloring(const Graph &graph, std::vector<std::size_t> &colors)

Compute vertex colors.

Parameters:
  • graph
  • colors

ParMETIS

C++ documentation for ParMETIS from dolfin/graph/ParMETIS.h:

class dolfin::ParMETIS

This class provides an interface to ParMETIS .

void dolfin::ParMETIS::adaptive_repartition(MPI_Comm mpi_comm, CSRGraph<T> &csr_graph, std::vector<int> &cell_partition)
Parameters:
  • mpi_comm
  • csr_graph
  • cell_partition
void dolfin::ParMETIS::compute_partition(const MPI_Comm mpi_comm, std::vector<int> &cell_partition, std::map<std::int64_t, std::vector<int>> &ghost_procs, const boost::multi_array<std::int64_t, 2> &cell_vertices, const std::size_t num_global_vertices, const CellType &cell_type, const std::string mode = "partition")

Compute cell partition from local mesh data. The output vector cell_partition contains the desired destination process numbers for each cell. Cells shared on multiple processes have an entry in ghost_procs pointing to the set of sharing process numbers. The mode argument determines which ParMETIS function is called. It can be one of “partition”, “adaptive_repartition” or “refine”. For meshes that have already been partitioned or are already well partitioned, it can be advantageous to use “adaptive_repartition” or “refine”.

Parameters:
  • mpi_comm
  • cell_partition
  • ghost_procs
  • cell_vertices
  • num_global_vertices
  • cell_type
  • mode
void dolfin::ParMETIS::partition(MPI_Comm mpi_comm, CSRGraph<T> &csr_graph, std::vector<int> &cell_partition, std::map<std::int64_t, std::vector<int>> &ghost_procs)
Parameters:
  • mpi_comm
  • csr_graph
  • cell_partition
  • ghost_procs
void dolfin::ParMETIS::refine(MPI_Comm mpi_comm, CSRGraph<T> &csr_graph, std::vector<int> &cell_partition)
Parameters:
  • mpi_comm
  • csr_graph
  • cell_partition

SCOTCH

C++ documentation for SCOTCH from dolfin/graph/SCOTCH.h:

class dolfin::SCOTCH

This class provides an interface to SCOTCH-PT (parallel version)

std::vector<int> dolfin::SCOTCH::compute_gps(const Graph &graph, std::size_t num_passes = 5)

Compute reordering (map[old] -> new) using Gibbs-Poole-Stockmeyer (GPS) re-ordering

Parameters:
  • graph – (Graph) Input graph
  • num_passes – (std::size_t) Number of passes to use in GPS algorithm
Returns:

std::vector<int> Mapping from old to new nodes

void dolfin::SCOTCH::compute_partition(const MPI_Comm mpi_comm, std::vector<int> &cell_partition, std::map<std::int64_t, std::vector<int>> &ghost_procs, const boost::multi_array<std::int64_t, 2> &cell_vertices, const std::vector<std::size_t> &cell_weight, const std::int64_t num_global_vertices, const std::int64_t num_global_cells, const CellType &cell_type)

Compute cell partition from local mesh data. The vector cell_partition contains the desired destination process numbers for each cell. Cells shared on multiple processes have an entry in ghost_procs pointing to the set of sharing process numbers.

Parameters:
  • mpi_comm – (MPI_Comm)
  • cell_partition – (std::vector<int>)
  • ghost_procs – (std::map<std::int64_t, std::vector<int>>)
  • cell_vertices – (const boost::multi_array<std::int64_t, 2>)
  • cell_weight – (const std::vector<std::size_t>)
  • num_global_vertices – (const std::int64_t)
  • num_global_cells – (const std::int64_t)
  • cell_type – (const CellType )
std::vector<int> dolfin::SCOTCH::compute_reordering(const Graph &graph, std::string scotch_strategy = "")

Compute graph re-ordering

Parameters:
  • graph – (Graph) Input graph
  • scotch_strategy – (string) SCOTCH parameters
Returns:

std::vector<int> Mapping from old to new nodes

void dolfin::SCOTCH::compute_reordering(const Graph &graph, std::vector<int> &permutation, std::vector<int> &inverse_permutation, std::string scotch_strategy = "")

Compute graph re-ordering

Parameters:
  • graph – (Graph)
  • permutation – (std::vector<int>)
  • inverse_permutation – (std::vector<int>)
  • scotch_strategy – (std::string)
void dolfin::SCOTCH::partition(const MPI_Comm mpi_comm, CSRGraph<T> &local_graph, const std::vector<std::size_t> &node_weights, const std::set<std::int64_t> &ghost_vertices, const std::size_t num_global_vertices, std::vector<int> &cell_partition, std::map<std::int64_t, std::vector<int>> &ghost_procs)
Parameters:
  • mpi_comm
  • local_graph
  • node_weights
  • ghost_vertices
  • num_global_vertices
  • cell_partition
  • ghost_procs

ZoltanInterface

C++ documentation for ZoltanInterface from dolfin/graph/ZoltanInterface.h:

class dolfin::ZoltanInterface

This class colors a graph using Zoltan (part of Trilinos). It is designed to work on a single process.

std::size_t dolfin::ZoltanInterface::compute_local_vertex_coloring(const Graph &graph, std::vector<std::size_t> &colors)

Compute vertex colors.

Parameters:
  • graph
  • colors

ZoltanGraphInterface

C++ documentation for ZoltanGraphInterface from dolfin/graph/ZoltanInterface.h:

class dolfin::ZoltanInterface::ZoltanGraphInterface
dolfin::ZoltanInterface::ZoltanGraphInterface::ZoltanGraphInterface(const Graph &graph)

Constructor.

Parameters:graph
void dolfin::ZoltanInterface::ZoltanGraphInterface::get_all_edges(void *data, int num_gid_entries, int num_lid_entries, int num_obj, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int *num_edges, ZOLTAN_ID_PTR nbor_global_id, int *nbor_procs, int wgt_dim, float *ewgts, int *ierr)
Parameters:
  • data
  • num_gid_entries
  • num_lid_entries
  • num_obj
  • global_ids
  • local_ids
  • num_edges
  • nbor_global_id
  • nbor_procs
  • wgt_dim
  • ewgts
  • ierr
void dolfin::ZoltanInterface::ZoltanGraphInterface::get_number_edges(void *data, int num_gid_entries, int num_lid_entries, int num_obj, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int *num_edges, int *ierr)
Parameters:
  • data
  • num_gid_entries
  • num_lid_entries
  • num_obj
  • global_ids
  • local_ids
  • num_edges
  • ierr
int dolfin::ZoltanInterface::ZoltanGraphInterface::get_number_of_objects(void *data, int *ierr)
Parameters:
  • data
  • ierr
void dolfin::ZoltanInterface::ZoltanGraphInterface::get_object_list(void *data, int sizeGID, int sizeLID, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int wgt_dim, float *obj_wgts, int *ierr)
Parameters:
  • data
  • sizeGID
  • sizeLID
  • globalID
  • localID
  • wgt_dim
  • obj_wgts
  • ierr
void dolfin::ZoltanInterface::ZoltanGraphInterface::num_vertex_edges(unsigned int *num_edges) const

Number of edges from each vertex.

Parameters:num_edges