Go to the documentation of this file.
37 #include <boost/config.hpp>
40 #include <boost/graph/adjacency_list.hpp>
41 #include <boost/graph/sloan_ordering.hpp>
42 #include <boost/graph/properties.hpp>
43 #include <boost/graph/bandwidth.hpp>
44 #include <boost/graph/profile.hpp>
45 #include <boost/graph/wavefront.hpp>
49 using namespace boost;
53 typedef adjacency_list
75 typedef graph_traits<Graph>::vertex_descriptor
Vertex;
76 typedef graph_traits<Graph>::vertices_size_type
size_type;
93 Foam::SloanRenumber::SloanRenumber(
const dictionary& renumberDict)
98 renumberDict.optionalSubDict
101 ).getOrDefault(
"reverse", false)
120 if (pp.coupled() && !isA<processorPolyPatch>(pp))
141 && !isA<processorPolyPatch>(pp)
142 && refCast<const coupledPolyPatch>(pp).owner()
145 label bFacei = pp.offset();
147 for (
const label ownCelli : pp.faceCells())
149 const label nbrCelli = nbr[bFacei];
152 if (ownCelli < nbrCelli)
154 add_edge(ownCelli, nbrCelli,
G);
158 add_edge(nbrCelli, ownCelli,
G);
166 graph_traits<Graph>::vertex_iterator ui, ui_end;
170 for (boost::tie(ui, ui_end) =
vertices(
G); ui != ui_end; ++ui)
171 deg[*ui] = degree(*ui,
G);
177 std::vector<Vertex> sloan_order(num_vertices(
G));
183 get(vertex_color,
G),
185 get(vertex_priority,
G)
188 labelList orderedToOld(sloan_order.size());
191 orderedToOld[
c] = index_map[sloan_order[
c]];
209 Graph G(cellCells.size());
213 const labelList& nbrs = cellCells[celli];
218 add_edge(celli, nbrs[i],
G);
224 graph_traits<Graph>::vertex_iterator ui, ui_end;
228 for (boost::tie(ui, ui_end) =
vertices(
G); ui != ui_end; ++ui)
229 deg[*ui] = degree(*ui,
G);
235 std::vector<Vertex> sloan_order(num_vertices(
G));
241 get(vertex_color,
G),
243 get(vertex_priority,
G)
246 labelList orderedToOld(sloan_order.size());
249 orderedToOld[
c] = index_map[sloan_order[
c]];
void reverse(UList< T > &list, const label n)
const dimensionedScalar G
Newtonian constant of gravitation.
A List obtained as a section of another List.
graph_traits< Graph >::vertex_descriptor Vertex
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited, i.e.
Mesh consisting of general polyhedral cells.
#define forAll(list, i)
Loop across all elements in list.
adjacency_list< setS, vecS, undirectedS, property< vertex_color_t, default_color_type, property< vertex_degree_t, Foam::label, property< vertex_priority_t, Foam::scalar > > >> Graph
label nCells() const noexcept
Number of mesh cells.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A patch is a list of labels that address the faces in the global face list.
virtual const labelList & faceOwner() const
Return face owner.
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
graph_traits< Graph >::vertices_size_type size_type
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
pointField vertices(const blockVertexList &bvl)
Macros for easy insertion into run-time selection tables.
Abstract base class for renumbering.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const dimensionedScalar c
Speed of light in a vacuum.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
virtual const labelList & faceNeighbour() const
Return face neighbour.