37#pragma clang diagnostic ignored "-Wdeprecated-builtins"
38#pragma clang diagnostic ignored "-Wdeprecated-declarations"
40#include <boost/config.hpp>
43#include <boost/graph/adjacency_list.hpp>
44#include <boost/graph/sloan_ordering.hpp>
45#include <boost/graph/properties.hpp>
46#include <boost/graph/bandwidth.hpp>
47#include <boost/graph/profile.hpp>
48#include <boost/graph/wavefront.hpp>
77typedef graph_traits<Graph>::vertex_descriptor
Vertex;
78typedef graph_traits<Graph>::vertices_size_type
size_type;
111 .getOrDefault(
"reverse", false)
123 using namespace Foam;
126 graph_traits<Graph>::vertex_iterator ui, ui_end;
129 property_map<Graph,vertex_degree_t>::type deg = get(vertex_degree, G);
130 for (boost::tie(ui, ui_end) =
vertices(G); ui != ui_end; ++ui)
132 deg[*ui] = degree(*ui, G);
136 property_map<Graph, vertex_index_t>::type index_map =
get(vertex_index, G);
139 std::vector<Vertex> sloan_order(num_vertices(G));
145 get(vertex_color, G),
147 get(vertex_priority, G)
150 labelList orderedToOld(sloan_order.size());
153 orderedToOld[
c] = index_map[sloan_order[
c]];
206 label bFacei =
pp.offset();
208 for (
const label ownCelli :
pp.faceCells())
210 const label nbrCelli = nbr[bFacei];
213 if (ownCelli < nbrCelli)
215 add_edge(ownCelli, nbrCelli, G);
219 add_edge(nbrCelli, ownCelli, G);
225 return renumberImpl(G, reverse_);
238 const auto& neighbours = cellCells[celli];
240 for (
const label nbr : neighbours)
244 add_edge(celli, nbr, G);
249 return renumberImpl(G, reverse_);
262 const auto& neighbours = cellCells[celli];
264 for (
const label nbr : neighbours)
268 add_edge(celli, nbr, G);
273 return renumberImpl(G, reverse_);
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
graph_traits< Graph >::vertices_size_type size_type
graph_traits< Graph >::vertex_descriptor Vertex
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A packed storage of objects of type <T> using an offset table for access.
label size() const noexcept
The primary size (the number of rows/sublists).
Sloan renumbering algorithm.
void reverse(bool on) noexcept
Toggle reverse on/off.
virtual labelList renumber(const polyMesh &mesh) const
Return the cell visit order (from ordered back to original cell id) using the mesh to determine the c...
SloanRenumber(const bool reverse=false)
Default construct, optionally with reverse.
A non-owning sub-view of a List (allocated or unallocated storage).
void size(const label n)
Older name for setAddressableSize.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Mesh consisting of general polyhedral cells.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
virtual const labelList & faceOwner() const
Return face owner.
virtual const labelList & faceNeighbour() const
Return face neighbour.
A patch is a list of labels that address the faces in the global face list.
label nInternalFaces() const noexcept
Number of internal faces.
label nCells() const noexcept
Number of mesh cells.
Abstract base class for renumbering.
renumberMethod()
Default construct.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const dimensionedScalar c
Speed of light in a vacuum.
const dimensionedScalar G
Newtonian constant of gravitation.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
pointField vertices(const blockVertexList &bvl)
#define forAll(list, i)
Loop across all elements in list.