SloanRenumber.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2017 OpenFOAM Foundation
9  Copyright (C) 2020-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27  Adapted from Boost graph/example/sloan_ordering.cpp
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #include "SloanRenumber.H"
33 #include "processorPolyPatch.H"
34 #include "syncTools.H"
35 
36 #include <boost/config.hpp>
37 #include <vector>
38 #include <iostream>
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>
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 using namespace boost;
49 
50 //Defining the graph type
51 typedef adjacency_list
52 <
53  setS,
54  vecS,
55  undirectedS,
56  property
57  <
58  vertex_color_t,
59  default_color_type,
60  property
61  <
62  vertex_degree_t,
63  Foam::label,
64  property
65  <
66  vertex_priority_t,
67  Foam::scalar
68  >
69  >
70  >
71 > Graph;
72 
73 typedef graph_traits<Graph>::vertex_descriptor Vertex;
74 typedef graph_traits<Graph>::vertices_size_type size_type;
75 
76 
77 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
78 
79 namespace Foam
80 {
81  defineTypeNameAndDebug(SloanRenumber, 0);
82 
84  (
85  renumberMethod,
86  SloanRenumber,
88  );
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
93 
95 :
97  reverse_(reverse)
98 {}
99 
100 
102 :
104  reverse_
105  (
106  dict.optionalSubDict(typeName + "Coeffs")
107  .getOrDefault("reverse", false)
108  )
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
113 
114 namespace
115 {
116 
117 Foam::labelList renumberImpl(Graph& G, const bool useReverse)
118 {
119  using namespace Foam;
120 
121  //Creating two iterators over the vertices
122  graph_traits<Graph>::vertex_iterator ui, ui_end;
123 
124  //Creating a property_map with the degrees of the degrees of each vertex
125  property_map<Graph,vertex_degree_t>::type deg = get(vertex_degree, G);
126  for (boost::tie(ui, ui_end) = vertices(G); ui != ui_end; ++ui)
127  {
128  deg[*ui] = degree(*ui, G);
129  }
130 
131  //Creating a property_map for the indices of a vertex
132  property_map<Graph, vertex_index_t>::type index_map = get(vertex_index, G);
133 
134  //Creating a vector of vertices
135  std::vector<Vertex> sloan_order(num_vertices(G));
136 
137  sloan_ordering
138  (
139  G,
140  sloan_order.begin(),
141  get(vertex_color, G),
142  make_degree_map(G),
143  get(vertex_priority, G)
144  );
145 
146  labelList orderedToOld(sloan_order.size());
147  forAll(orderedToOld, c)
148  {
149  orderedToOld[c] = index_map[sloan_order[c]];
150  }
151 
152  if (useReverse)
153  {
154  Foam::reverse(orderedToOld);
155  }
156 
157  return orderedToOld;
158 }
160 } // End anonymous namespace
161 
162 
163 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164 
166 (
167  const polyMesh& mesh
168 ) const
169 {
170  // Construct graph : faceOwner + connections across cyclics.
171 
172  // Determine neighbour cell
173  labelList nbr(mesh.nBoundaryFaces(), -1);
174  for (const polyPatch& pp : mesh.boundaryMesh())
175  {
176  if (pp.coupled() && !isA<processorPolyPatch>(pp))
177  {
178  SubList<label>(nbr, pp.size(), pp.offset()) = pp.faceCells();
179  }
180  }
182 
183 
184  Graph G(mesh.nCells());
185 
186  // Add internal faces
187  for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
188  {
189  add_edge(mesh.faceOwner()[facei], mesh.faceNeighbour()[facei], G);
190  }
191 
192  // Add cyclics
193  for (const polyPatch& pp : mesh.boundaryMesh())
194  {
195  if
196  (
197  pp.coupled()
198  && !isA<processorPolyPatch>(pp)
199  && refCast<const coupledPolyPatch>(pp).owner()
200  )
201  {
202  label bFacei = pp.offset();
203 
204  for (const label ownCelli : pp.faceCells())
205  {
206  const label nbrCelli = nbr[bFacei];
207  ++bFacei;
208 
209  if (ownCelli < nbrCelli)
210  {
211  add_edge(ownCelli, nbrCelli, G);
212  }
213  else
214  {
215  add_edge(nbrCelli, ownCelli, G);
216  }
217  }
218  }
219  }
220 
221  return renumberImpl(G, reverse_);
222 }
223 
224 
226 (
227  const CompactListList<label>& cellCells
228 ) const
229 {
230  Graph G(cellCells.size());
231 
232  forAll(cellCells, celli)
233  {
234  const auto& neighbours = cellCells[celli];
235 
236  for (const label nbr : neighbours)
237  {
238  if (celli < nbr)
239  {
240  add_edge(celli, nbr, G);
241  }
242  }
243  }
244 
245  return renumberImpl(G, reverse_);
246 }
247 
248 
250 (
251  const labelListList& cellCells
252 ) const
253 {
254  Graph G(cellCells.size());
255 
256  forAll(cellCells, celli)
257  {
258  const auto& neighbours = cellCells[celli];
259 
260  for (const label nbr : neighbours)
261  {
262  if (celli < nbr)
263  {
264  add_edge(celli, nbr, G);
265  }
266  }
267  }
268 
269  return renumberImpl(G, reverse_);
270 }
271 
272 
273 // ************************************************************************* //
graph_traits< Graph >::vertex_descriptor Vertex
Definition: SloanRenumber.C:66
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
type
Types of root.
Definition: Roots.H:52
label size() const noexcept
The primary size (the number of rows/sublists)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
const dimensionedScalar G
Newtonian constant of gravitation.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const bool parRun=UPstream::parRun())
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:524
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.
Definition: SloanRenumber.C:87
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
pointField vertices(const blockVertexList &bvl)
Abstract base class for renumbering.
A List obtained as a section of another List.
Definition: SubList.H:50
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:142
dynamicFvMesh & mesh
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
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
Definition: SloanRenumber.C:64
label nInternalFaces() const noexcept
Number of internal faces.
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition: UListI.H:521
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:67
A packed storage of objects of type <T> using an offset table for access.
label nCells() const noexcept
Number of mesh cells.
const dimensionedScalar c
Speed of light in a vacuum.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.