Go to the documentation of this file.
34 #include <boost/config.hpp>
37 #include <boost/graph/adjacency_list.hpp>
38 #include <boost/graph/sloan_ordering.hpp>
39 #include <boost/graph/properties.hpp>
40 #include <boost/graph/bandwidth.hpp>
41 #include <boost/graph/profile.hpp>
42 #include <boost/graph/wavefront.hpp>
46 using namespace boost;
50 typedef adjacency_list
72 typedef graph_traits<Graph>::vertex_descriptor
Vertex;
73 typedef graph_traits<Graph>::vertices_size_type
size_type;
95 renumberDict.
found(typeName +
"Coeffs")
96 ?
Switch(renumberDict.subDict(typeName +
"Coeffs").
lookup(
"reverse"))
118 if (pbm[patchI].coupled() && !isA<processorPolyPatch>(pbm[patchI]))
125 ).
assign(pbm[patchI].faceCells());
143 pbm[patchI].coupled()
144 && !isA<processorPolyPatch>(pbm[patchI])
145 && refCast<const coupledPolyPatch>(pbm[patchI]).owner()
148 const labelUList& faceCells = pbm[patchI].faceCells();
152 label nbrCellI = nbr[bFaceI];
154 if (faceCells[i] < nbrCellI)
156 add_edge(faceCells[i], nbrCellI,
G);
160 add_edge(nbrCellI, faceCells[i],
G);
168 graph_traits<Graph>::vertex_iterator ui, ui_end;
172 for (boost::tie(ui, ui_end) = vertices(
G); ui != ui_end; ++ui)
173 deg[*ui] = degree(*ui,
G);
179 std::vector<Vertex> sloan_order(num_vertices(
G));
185 get(vertex_color,
G),
187 get(vertex_priority,
G)
190 labelList orderedToOld(sloan_order.size());
193 orderedToOld[
c] = index_map[sloan_order[
c]];
215 const labelList& nbrs = cellCells[cellI];
220 add_edge(cellI, nbrs[i],
G);
226 graph_traits<Graph>::vertex_iterator ui, ui_end;
230 for (boost::tie(ui, ui_end) = vertices(
G); ui != ui_end; ++ui)
231 deg[*ui] = degree(*ui,
G);
237 std::vector<Vertex> sloan_order(num_vertices(
G));
243 get(vertex_color,
G),
245 get(vertex_priority,
G)
248 labelList orderedToOld(sloan_order.size());
251 orderedToOld[
c] = index_map[sloan_order[
c]];
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
#define forAll(list, i)
Loop across all elements in list.
const dimensionedScalar G
Newtonian constant of gravitation.
SloanRenumber(const SloanRenumber &)
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.
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
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
virtual const labelList & faceOwner() const
Return face owner.
void assign(const UList< T > &)
Assign elements to those from UList.
graph_traits< Graph >::vertices_size_type size_type
A list of keyword definitions, which are a keyword followed by any number of values (e....
label nInternalFaces() const
Macros for easy insertion into run-time selection tables.
defineTypeNameAndDebug(SloanRenumber, 0)
Abstract base class for renumbering.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
const dimensionedScalar c
Speed of light in a vacuum.
void size(const label)
Override size to be inconsistent with allocated storage.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
virtual const labelList & faceNeighbour() const
Return face neighbour.
stressControl lookup("compactNormalStress") >> compactNormalStress
void reverse(UList< T > &, const label n)