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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26  Adapted from Boost graph/example/sloan_ordering.cpp
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "SloanRenumber.H"
32 #include "decompositionMethod.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 using namespace std;
50 
51 //Defining the graph type
52 typedef adjacency_list
53 <
54  setS,
55  vecS,
56  undirectedS,
57  property
58  <
59  vertex_color_t,
60  default_color_type,
61  property
62  <
63  vertex_degree_t,
65  property
66  <
67  vertex_priority_t,
68  Foam::scalar
69  >
70  >
71  >
73 
74 typedef graph_traits<Graph>::vertex_descriptor Vertex;
75 typedef graph_traits<Graph>::vertices_size_type size_type;
76 
77 namespace Foam
78 {
79  defineTypeNameAndDebug(SloanRenumber, 0);
80 
82  (
83  renumberMethod,
84  SloanRenumber,
85  dictionary
86  );
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91 
92 Foam::SloanRenumber::SloanRenumber(const dictionary& renumberDict)
93 :
94  renumberMethod(renumberDict),
95  reverse_
96  (
97  renumberDict.optionalSubDict
98  (
99  typeName + "Coeffs"
100  ).lookupOrDefault("reverse", false)
101  )
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 (
109  const polyMesh& mesh,
110  const pointField& points
111 ) const
112 {
113  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
114 
115  // Construct graph : faceOwner + connections across cyclics.
116 
117  // Determine neighbour cell
118  labelList nbr(mesh.nBoundaryFaces(), -1);
119  forAll(pbm, patchi)
120  {
121  if (pbm[patchi].coupled() && !isA<processorPolyPatch>(pbm[patchi]))
122  {
124  (
125  nbr,
126  pbm[patchi].size(),
127  pbm[patchi].start()-mesh.nInternalFaces()
128  ) = pbm[patchi].faceCells();
129  }
130  }
132 
133 
134  Graph G(mesh.nCells());
135 
136  // Add internal faces
137  forAll(mesh.faceNeighbour(), facei)
138  {
139  add_edge(mesh.faceOwner()[facei], mesh.faceNeighbour()[facei], G);
140  }
141  // Add cyclics
142  forAll(pbm, patchi)
143  {
144  if
145  (
146  pbm[patchi].coupled()
147  && !isA<processorPolyPatch>(pbm[patchi])
148  && refCast<const coupledPolyPatch>(pbm[patchi]).owner()
149  )
150  {
151  const labelUList& faceCells = pbm[patchi].faceCells();
152  forAll(faceCells, i)
153  {
154  label bFacei = pbm[patchi].start()+i-mesh.nInternalFaces();
155  label nbrCelli = nbr[bFacei];
156 
157  if (faceCells[i] < nbrCelli)
158  {
159  add_edge(faceCells[i], nbrCelli, G);
160  }
161  else
162  {
163  add_edge(nbrCelli, faceCells[i], G);
164  }
165  }
166  }
167  }
168 
169 
170  //Creating two iterators over the vertices
171  graph_traits<Graph>::vertex_iterator ui, ui_end;
172 
173  //Creating a property_map with the degrees of the degrees of each vertex
174  property_map<Graph,vertex_degree_t>::type deg = get(vertex_degree, G);
175  for (boost::tie(ui, ui_end) = vertices(G); ui != ui_end; ++ui)
176  deg[*ui] = degree(*ui, G);
177 
178  //Creating a property_map for the indices of a vertex
179  property_map<Graph, vertex_index_t>::type index_map = get(vertex_index, G);
180 
181  //Creating a vector of vertices
182  std::vector<Vertex> sloan_order(num_vertices(G));
183 
184  sloan_ordering
185  (
186  G,
187  sloan_order.begin(),
188  get(vertex_color, G),
189  make_degree_map(G),
190  get(vertex_priority, G)
191  );
192 
193  labelList orderedToOld(sloan_order.size());
194  forAll(orderedToOld, c)
195  {
196  orderedToOld[c] = index_map[sloan_order[c]];
197  }
198 
199  if (reverse_)
200  {
201  reverse(orderedToOld);
202  }
203 
204  return orderedToOld;
205 }
206 
207 
209 (
210  const labelListList& cellCells,
211  const pointField& points
212 ) const
213 {
214  Graph G(cellCells.size());
215 
216  forAll(cellCells, celli)
217  {
218  const labelList& nbrs = cellCells[celli];
219  forAll(nbrs, i)
220  {
221  if (nbrs[i] > celli)
222  {
223  add_edge(celli, nbrs[i], G);
224  }
225  }
226  }
227 
228  //Creating two iterators over the vertices
229  graph_traits<Graph>::vertex_iterator ui, ui_end;
230 
231  //Creating a property_map with the degrees of the degrees of each vertex
232  property_map<Graph,vertex_degree_t>::type deg = get(vertex_degree, G);
233  for (boost::tie(ui, ui_end) = vertices(G); ui != ui_end; ++ui)
234  deg[*ui] = degree(*ui, G);
235 
236  //Creating a property_map for the indices of a vertex
237  property_map<Graph, vertex_index_t>::type index_map = get(vertex_index, G);
238 
239  //Creating a vector of vertices
240  std::vector<Vertex> sloan_order(num_vertices(G));
241 
242  sloan_ordering
243  (
244  G,
245  sloan_order.begin(),
246  get(vertex_color, G),
247  make_degree_map(G),
248  get(vertex_priority, G)
249  );
250 
251  labelList orderedToOld(sloan_order.size());
252  forAll(orderedToOld, c)
253  {
254  orderedToOld[c] = index_map[sloan_order[c]];
255  }
256 
257  if (reverse_)
258  {
259  reverse(orderedToOld);
260  }
261 
262  return orderedToOld;
263 }
264 
265 
266 // ************************************************************************* //
Foam::reverse
void reverse(UList< T > &list, const label n)
Definition: UListI.H:396
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Vertex
graph_traits< Graph >::vertex_descriptor Vertex
Definition: SloanRenumber.C:74
SloanRenumber.H
Foam::syncTools::swapBoundaryFaceList
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:439
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
decompositionMethod.H
syncTools.H
Foam::polyBoundaryMesh::start
label start() const
The start label of the boundary faces in the polyMesh face list.
Definition: polyBoundaryMesh.C:638
Foam::SloanRenumber::renumber
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited, i.e.
Definition: SloanRenumber.H:90
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Graph
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:72
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::label
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:62
Foam::Field< vector >
addToRunTimeSelectionTable
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to hash-table of functions with typename as the key.
Definition: addToRunTimeSelectionTable.H:37
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1076
size_type
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:75
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
processorPolyPatch.H
Foam::primitiveMesh::nBoundaryFaces
label nBoundaryFaces() const
Number of boundary faces (== nFaces - nInternalFaces)
Definition: primitiveMeshI.H:84
Foam::vertices
pointField vertices(const blockVertexList &bvl)
Definition: blockVertexList.H:49
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::renumberMethod
Abstract base class for renumbering.
Definition: renumberMethod.H:50
Foam::List< label >
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
Foam::UList< label >
defineTypeNameAndDebug
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1082