19#ifndef OPM_GRAPHCOLORING_HEADER_INCLUDED
20#define OPM_GRAPHCOLORING_HEADER_INCLUDED
37std::size_t colorGraphWelshPowell(
const Graph& graph,
38 std::deque<typename Graph::VertexDescriptor>& orderedVertices,
39 std::vector<int>& colors,
40 int color,
int noVertices)
42 std::vector<int> forbidden(noVertices,
false);
43 std::size_t noColored = 0;
45 for(
auto vertex = orderedVertices.begin(),
46 vertexEnd = orderedVertices.end();
47 vertex != vertexEnd; ++vertex)
50 while(vertex != vertexEnd && forbidden[*vertex])
52 if ( vertex == vertexEnd )
58 colors[*vertex] = color;
61 for(
auto edge = graph.beginEdges(*vertex), endEdge = graph.endEdges(*vertex);
62 edge != endEdge; ++edge)
64 forbidden[edge.target()] =
true;
68 using Vertex =
typename Graph::VertexDescriptor;
69 auto newEnd = std::remove_if(orderedVertices.begin(), orderedVertices.end(),
70 [&forbidden](
const Vertex& vertex)
72 return !forbidden[vertex];
74 orderedVertices.resize(newEnd-orderedVertices.begin());
77template<
class Graph,
class Functor>
78std::size_t breadthFirstSearch(
const Graph& graph,
typename Graph::VertexDescriptor root,
81 std::vector<int> visited(graph.maxVertex() + 1,
false);
82 using Vertex =
typename Graph::VertexDescriptor;
83 std::queue<Vertex> nextVertices;
84 std::size_t noVisited = 0;
85 nextVertices.push(root);
88 while( !nextVertices.empty() )
90 auto current = nextVertices.front();
91 for(
auto edge = graph.beginEdges(current),
92 endEdge = graph.endEdges(current);
93 edge != endEdge; ++edge)
95 if ( ! visited[edge.target()] )
97 visited[edge.target()] =
true;
98 nextVertices.push(edge.target());
99 functor(edge.target());
117std::tuple<std::vector<int>, int, std::vector<std::size_t> >
120 using Vertex =
typename Graph::VertexDescriptor;
163 return degrees[v1] > degrees[v2];
185std::vector<std::size_t>
191 std::vector<std::size_t>
indices(
graph.maxVertex() + 1);
205std::vector<std::size_t>
209 typename Graph::VertexDescriptor
root)
212 const auto notVisitedTag = std::numeric_limits<std::size_t>::max();
214 using Vertex =
typename Graph::VertexDescriptor;
Definition AquiferInterface.hpp:35
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27
std::vector< std::size_t > reorderVerticesPreserving(const std::vector< int > &colors, int noColors, const std::vector< std::size_t > &verticesPerColor, const Graph &graph)
! Reorder colored graph preserving order of vertices with the same color.
Definition GraphColoring.hpp:186
std::vector< std::size_t > reorderVerticesSpheres(const std::vector< int > &colors, int noColors, const std::vector< std::size_t > &verticesPerColor, const Graph &graph, typename Graph::VertexDescriptor root)
! Reorder Vetrices in spheres
Definition GraphColoring.hpp:206
std::tuple< std::vector< int >, int, std::vector< std::size_t > > colorVerticesWelshPowell(const Graph &graph)
Color the vertices of graph.
Definition GraphColoring.hpp:118