Go to the documentation of this file.
36 #include <boost/config.hpp>
39 #include <boost/graph/adjacency_list.hpp>
40 #include <boost/graph/sloan_ordering.hpp>
41 #include <boost/graph/properties.hpp>
42 #include <boost/graph/bandwidth.hpp>
43 #include <boost/graph/profile.hpp>
44 #include <boost/graph/wavefront.hpp>
48 using namespace boost;
52 typedef adjacency_list
74 typedef graph_traits<Graph>::vertex_descriptor
Vertex;
75 typedef graph_traits<Graph>::vertices_size_type
size_type;
92 Foam::SloanRenumber::SloanRenumber(
const dictionary& renumberDict)
97 renumberDict.optionalSubDict
100 ).lookupOrDefault(
"reverse", false)
121 if (pbm[patchi].coupled() && !isA<processorPolyPatch>(pbm[patchi]))
146 pbm[patchi].coupled()
147 && !isA<processorPolyPatch>(pbm[patchi])
148 && refCast<const coupledPolyPatch>(pbm[patchi]).owner()
155 label nbrCelli = nbr[bFacei];
171 graph_traits<Graph>::vertex_iterator ui, ui_end;
175 for (boost::tie(ui, ui_end) =
vertices(
G); ui != ui_end; ++ui)
176 deg[*ui] = degree(*ui,
G);
182 std::vector<Vertex> sloan_order(num_vertices(
G));
188 get(vertex_color,
G),
190 get(vertex_priority,
G)
193 labelList orderedToOld(sloan_order.size());
196 orderedToOld[
c] = index_map[sloan_order[
c]];
214 Graph G(cellCells.size());
218 const labelList& nbrs = cellCells[celli];
223 add_edge(celli, nbrs[i],
G);
229 graph_traits<Graph>::vertex_iterator ui, ui_end;
233 for (boost::tie(ui, ui_end) =
vertices(
G); ui != ui_end; ++ui)
234 deg[*ui] = degree(*ui,
G);
240 std::vector<Vertex> sloan_order(num_vertices(
G));
246 get(vertex_color,
G),
248 get(vertex_priority,
G)
251 labelList orderedToOld(sloan_order.size());
254 orderedToOld[
c] = index_map[sloan_order[
c]];
void reverse(UList< T > &list, const label n)
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
const dimensionedScalar G
Newtonian constant of gravitation.
label nInternalFaces() const
Number of internal faces.
A List obtained as a section of another List.
graph_traits< Graph >::vertex_descriptor Vertex
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
label start() const
The start label of the boundary faces in the polyMesh face list.
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
Number of mesh cells.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to hash-table of functions with typename as the key.
virtual const labelList & faceOwner() const
Return face owner.
graph_traits< Graph >::vertices_size_type size_type
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
label nBoundaryFaces() const
Number of boundary faces (== nFaces - nInternalFaces)
pointField vertices(const blockVertexList &bvl)
Macros for easy insertion into run-time selection tables.
Abstract base class for renumbering.
label ListType::const_reference const label start
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const dimensionedScalar c
Speed of light in a vacuum.
Smooth ATC in cells next to a set of patches supplied by type.
virtual const labelList & faceNeighbour() const
Return face neighbour.