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 | Copyright (C) 2012-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24  Adapted from Boost graph/example/sloan_ordering.cpp
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "SloanRenumber.H"
30 #include "decompositionMethod.H"
31 #include "processorPolyPatch.H"
32 #include "syncTools.H"
33 
34 #include <boost/config.hpp>
35 #include <vector>
36 #include <iostream>
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>
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 using namespace boost;
47 using namespace std;
48 
49 //Defining the graph type
50 typedef adjacency_list
51 <
52  setS,
53  vecS,
54  undirectedS,
55  property
56  <
57  vertex_color_t,
58  default_color_type,
59  property
60  <
61  vertex_degree_t,
63  property
64  <
65  vertex_priority_t,
66  Foam::scalar
67  >
68  >
69  >
71 
72 typedef graph_traits<Graph>::vertex_descriptor Vertex;
73 typedef graph_traits<Graph>::vertices_size_type size_type;
74 
75 namespace Foam
76 {
77  defineTypeNameAndDebug(SloanRenumber, 0);
78 
80  (
81  renumberMethod,
82  SloanRenumber,
83  dictionary
84  );
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89 
90 Foam::SloanRenumber::SloanRenumber(const dictionary& renumberDict)
91 :
92  renumberMethod(renumberDict),
93  reverse_
94  (
95  renumberDict.found(typeName + "Coeffs")
96  ? Switch(renumberDict.subDict(typeName + "Coeffs").lookup("reverse"))
97  : Switch(false)
98  )
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
105 (
106  const polyMesh& mesh,
107  const pointField& points
108 ) const
109 {
110  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
111 
112  // Construct graph : faceOwner + connections across cyclics.
113 
114  // Determine neighbour cell
115  labelList nbr(mesh.nFaces()-mesh.nInternalFaces(), -1);
116  forAll(pbm, patchi)
117  {
118  if (pbm[patchi].coupled() && !isA<processorPolyPatch>(pbm[patchi]))
119  {
121  (
122  nbr,
123  pbm[patchi].size(),
124  pbm[patchi].start()-mesh.nInternalFaces()
125  ) = pbm[patchi].faceCells();
126  }
127  }
129 
130 
131  Graph G(mesh.nCells());
132 
133  // Add internal faces
134  forAll(mesh.faceNeighbour(), facei)
135  {
136  add_edge(mesh.faceOwner()[facei], mesh.faceNeighbour()[facei], G);
137  }
138  // Add cyclics
139  forAll(pbm, patchi)
140  {
141  if
142  (
143  pbm[patchi].coupled()
144  && !isA<processorPolyPatch>(pbm[patchi])
145  && refCast<const coupledPolyPatch>(pbm[patchi]).owner()
146  )
147  {
148  const labelUList& faceCells = pbm[patchi].faceCells();
149  forAll(faceCells, i)
150  {
151  label bFacei = pbm[patchi].start()+i-mesh.nInternalFaces();
152  label nbrCelli = nbr[bFacei];
153 
154  if (faceCells[i] < nbrCelli)
155  {
156  add_edge(faceCells[i], nbrCelli, G);
157  }
158  else
159  {
160  add_edge(nbrCelli, faceCells[i], G);
161  }
162  }
163  }
164  }
165 
166 
167  //Creating two iterators over the vertices
168  graph_traits<Graph>::vertex_iterator ui, ui_end;
169 
170  //Creating a property_map with the degrees of the degrees of each vertex
171  property_map<Graph,vertex_degree_t>::type deg = get(vertex_degree, G);
172  for (boost::tie(ui, ui_end) = vertices(G); ui != ui_end; ++ui)
173  deg[*ui] = degree(*ui, G);
174 
175  //Creating a property_map for the indices of a vertex
176  property_map<Graph, vertex_index_t>::type index_map = get(vertex_index, G);
177 
178  //Creating a vector of vertices
179  std::vector<Vertex> sloan_order(num_vertices(G));
180 
181  sloan_ordering
182  (
183  G,
184  sloan_order.begin(),
185  get(vertex_color, G),
186  make_degree_map(G),
187  get(vertex_priority, G)
188  );
189 
190  labelList orderedToOld(sloan_order.size());
191  forAll(orderedToOld, c)
192  {
193  orderedToOld[c] = index_map[sloan_order[c]];
194  }
195 
196  if (reverse_)
197  {
198  reverse(orderedToOld);
199  }
200 
201  return orderedToOld;
202 }
203 
204 
206 (
207  const labelListList& cellCells,
208  const pointField& points
209 ) const
210 {
211  Graph G(cellCells.size());
212 
213  forAll(cellCells, celli)
214  {
215  const labelList& nbrs = cellCells[celli];
216  forAll(nbrs, i)
217  {
218  if (nbrs[i] > celli)
219  {
220  add_edge(celli, nbrs[i], G);
221  }
222  }
223  }
224 
225  //Creating two iterators over the vertices
226  graph_traits<Graph>::vertex_iterator ui, ui_end;
227 
228  //Creating a property_map with the degrees of the degrees of each vertex
229  property_map<Graph,vertex_degree_t>::type deg = get(vertex_degree, G);
230  for (boost::tie(ui, ui_end) = vertices(G); ui != ui_end; ++ui)
231  deg[*ui] = degree(*ui, G);
232 
233  //Creating a property_map for the indices of a vertex
234  property_map<Graph, vertex_index_t>::type index_map = get(vertex_index, G);
235 
236  //Creating a vector of vertices
237  std::vector<Vertex> sloan_order(num_vertices(G));
238 
239  sloan_ordering
240  (
241  G,
242  sloan_order.begin(),
243  get(vertex_color, G),
244  make_degree_map(G),
245  get(vertex_priority, G)
246  );
247 
248  labelList orderedToOld(sloan_order.size());
249  forAll(orderedToOld, c)
250  {
251  orderedToOld[c] = index_map[sloan_order[c]];
252  }
253 
254  if (reverse_)
255  {
256  reverse(orderedToOld);
257  }
258 
259  return orderedToOld;
260 }
261 
262 
263 // ************************************************************************* //
graph_traits< Graph >::vertex_descriptor Vertex
Definition: SloanRenumber.C:72
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
const dimensionedScalar G
Newtonian constant of gravitation.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Macros for easy insertion into run-time selection tables.
label nCells() const
Abstract base class for renumbering.
A List obtained as a section of another List.
Definition: SubList.H:53
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:119
stressControl lookup("compactNormalStress") >> compactNormalStress
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:70
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
Foam::polyBoundaryMesh.
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited, i.e.
Definition: SloanRenumber.H:87
void reverse(UList< T > &, const label n)
Definition: UListI.H:322
label nFaces() const
label patchi
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
const dimensionedScalar c
Speed of light in a vacuum.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:430
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
label nInternalFaces() const
bool found
Namespace for OpenFOAM.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29