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-2021 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 "decompositionMethod.H"
34 #include "processorPolyPatch.H"
35 #include "syncTools.H"
36 
37 #include <boost/config.hpp>
38 #include <vector>
39 #include <iostream>
40 #include <boost/graph/adjacency_list.hpp>
41 #include <boost/graph/sloan_ordering.hpp>
42 #include <boost/graph/properties.hpp>
43 #include <boost/graph/bandwidth.hpp>
44 #include <boost/graph/profile.hpp>
45 #include <boost/graph/wavefront.hpp>
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 using namespace boost;
50 using namespace std;
51 
52 //Defining the graph type
53 typedef adjacency_list
54 <
55  setS,
56  vecS,
57  undirectedS,
58  property
59  <
60  vertex_color_t,
61  default_color_type,
62  property
63  <
64  vertex_degree_t,
65  Foam::label,
66  property
67  <
68  vertex_priority_t,
69  Foam::scalar
70  >
71  >
72  >
74 
75 typedef graph_traits<Graph>::vertex_descriptor Vertex;
76 typedef graph_traits<Graph>::vertices_size_type size_type;
77 
78 namespace Foam
79 {
80  defineTypeNameAndDebug(SloanRenumber, 0);
81 
83  (
84  renumberMethod,
85  SloanRenumber,
86  dictionary
87  );
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
93 Foam::SloanRenumber::SloanRenumber(const dictionary& renumberDict)
94 :
95  renumberMethod(renumberDict),
96  reverse_
97  (
98  renumberDict.optionalSubDict
99  (
100  typeName + "Coeffs"
101  ).getOrDefault("reverse", false)
102  )
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
109 (
110  const polyMesh& mesh,
111  const pointField& points
112 ) const
113 {
114  // Construct graph : faceOwner + connections across cyclics.
115 
116  // Determine neighbour cell
117  labelList nbr(mesh.nBoundaryFaces(), -1);
118  for (const polyPatch& pp : mesh.boundaryMesh())
119  {
120  if (pp.coupled() && !isA<processorPolyPatch>(pp))
121  {
122  SubList<label>(nbr, pp.size(), pp.offset()) = pp.faceCells();
123  }
124  }
126 
127 
128  Graph G(mesh.nCells());
129 
130  // Add internal faces
131  forAll(mesh.faceNeighbour(), facei)
132  {
133  add_edge(mesh.faceOwner()[facei], mesh.faceNeighbour()[facei], G);
134  }
135  // Add cyclics
136  for (const polyPatch& pp : mesh.boundaryMesh())
137  {
138  if
139  (
140  pp.coupled()
141  && !isA<processorPolyPatch>(pp)
142  && refCast<const coupledPolyPatch>(pp).owner()
143  )
144  {
145  label bFacei = pp.offset();
146 
147  for (const label ownCelli : pp.faceCells())
148  {
149  const label nbrCelli = nbr[bFacei];
150  ++bFacei;
151 
152  if (ownCelli < nbrCelli)
153  {
154  add_edge(ownCelli, nbrCelli, G);
155  }
156  else
157  {
158  add_edge(nbrCelli, ownCelli, G);
159  }
160  }
161  }
162  }
163 
164 
165  //Creating two iterators over the vertices
166  graph_traits<Graph>::vertex_iterator ui, ui_end;
167 
168  //Creating a property_map with the degrees of the degrees of each vertex
169  property_map<Graph,vertex_degree_t>::type deg = get(vertex_degree, G);
170  for (boost::tie(ui, ui_end) = vertices(G); ui != ui_end; ++ui)
171  deg[*ui] = degree(*ui, G);
172 
173  //Creating a property_map for the indices of a vertex
174  property_map<Graph, vertex_index_t>::type index_map = get(vertex_index, G);
175 
176  //Creating a vector of vertices
177  std::vector<Vertex> sloan_order(num_vertices(G));
178 
179  sloan_ordering
180  (
181  G,
182  sloan_order.begin(),
183  get(vertex_color, G),
184  make_degree_map(G),
185  get(vertex_priority, G)
186  );
187 
188  labelList orderedToOld(sloan_order.size());
189  forAll(orderedToOld, c)
190  {
191  orderedToOld[c] = index_map[sloan_order[c]];
192  }
193 
194  if (reverse_)
195  {
196  reverse(orderedToOld);
197  }
198 
199  return orderedToOld;
200 }
201 
202 
204 (
205  const labelListList& cellCells,
206  const pointField& points
207 ) const
208 {
209  Graph G(cellCells.size());
210 
211  forAll(cellCells, celli)
212  {
213  const labelList& nbrs = cellCells[celli];
214  forAll(nbrs, i)
215  {
216  if (nbrs[i] > celli)
217  {
218  add_edge(celli, nbrs[i], G);
219  }
220  }
221  }
222 
223  //Creating two iterators over the vertices
224  graph_traits<Graph>::vertex_iterator ui, ui_end;
225 
226  //Creating a property_map with the degrees of the degrees of each vertex
227  property_map<Graph,vertex_degree_t>::type deg = get(vertex_degree, G);
228  for (boost::tie(ui, ui_end) = vertices(G); ui != ui_end; ++ui)
229  deg[*ui] = degree(*ui, G);
230 
231  //Creating a property_map for the indices of a vertex
232  property_map<Graph, vertex_index_t>::type index_map = get(vertex_index, G);
233 
234  //Creating a vector of vertices
235  std::vector<Vertex> sloan_order(num_vertices(G));
236 
237  sloan_ordering
238  (
239  G,
240  sloan_order.begin(),
241  get(vertex_color, G),
242  make_degree_map(G),
243  get(vertex_priority, G)
244  );
245 
246  labelList orderedToOld(sloan_order.size());
247  forAll(orderedToOld, c)
248  {
249  orderedToOld[c] = index_map[sloan_order[c]];
250  }
251 
252  if (reverse_)
253  {
254  reverse(orderedToOld);
255  }
256 
257  return orderedToOld;
258 }
259 
260 
261 // ************************************************************************* //
Foam::reverse
void reverse(UList< T > &list, const label n)
Definition: UListI.H:449
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:54
Vertex
graph_traits< Graph >::vertex_descriptor Vertex
Definition: SloanRenumber.C:75
SloanRenumber.H
Foam::syncTools::swapBoundaryFaceList
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:445
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
decompositionMethod.H
syncTools.H
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:296
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:73
Foam::primitiveMesh::nCells
label nCells() const noexcept
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::Field< vector >
addToRunTimeSelectionTable
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Definition: addToRunTimeSelectionTable.H:42
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::primitiveMesh::nBoundaryFaces
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
Definition: primitiveMeshI.H:84
size_type
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:76
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
processorPolyPatch.H
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 >
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::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113