Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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// With older boost...
37#pragma clang diagnostic ignored "-Wdeprecated-builtins"
38#pragma clang diagnostic ignored "-Wdeprecated-declarations"
39
40#include <boost/config.hpp>
41#include <vector>
42#include <iostream>
43#include <boost/graph/adjacency_list.hpp>
44#include <boost/graph/sloan_ordering.hpp>
45#include <boost/graph/properties.hpp>
46#include <boost/graph/bandwidth.hpp>
47#include <boost/graph/profile.hpp>
48#include <boost/graph/wavefront.hpp>
49
50// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51
52using namespace boost;
53
54//Defining the graph type
55typedef adjacency_list
56<
57 setS,
58 vecS,
59 undirectedS,
60 property
61 <
62 vertex_color_t,
63 default_color_type,
64 property
65 <
66 vertex_degree_t,
67 Foam::label,
68 property
69 <
70 vertex_priority_t,
71 Foam::scalar
72 >
73 >
74 >
75> Graph;
76
77typedef graph_traits<Graph>::vertex_descriptor Vertex;
78typedef graph_traits<Graph>::vertices_size_type size_type;
79
81// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
82
83namespace Foam
84{
86
88 (
92 );
93}
94
95
96// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97
108 reverse_
109 (
110 dict.optionalSubDict(typeName + "Coeffs")
111 .getOrDefault("reverse", false)
112 )
113{}
114
115
116// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
117
118namespace
119{
120
121Foam::labelList renumberImpl(Graph& G, const bool useReverse)
122{
123 using namespace Foam;
124
125 //Creating two iterators over the vertices
126 graph_traits<Graph>::vertex_iterator ui, ui_end;
127
128 //Creating a property_map with the degrees of the degrees of each vertex
129 property_map<Graph,vertex_degree_t>::type deg = get(vertex_degree, G);
130 for (boost::tie(ui, ui_end) = vertices(G); ui != ui_end; ++ui)
131 {
132 deg[*ui] = degree(*ui, G);
133 }
134
135 //Creating a property_map for the indices of a vertex
136 property_map<Graph, vertex_index_t>::type index_map = get(vertex_index, G);
137
138 //Creating a vector of vertices
139 std::vector<Vertex> sloan_order(num_vertices(G));
140
141 sloan_ordering
142 (
143 G,
144 sloan_order.begin(),
145 get(vertex_color, G),
146 make_degree_map(G),
147 get(vertex_priority, G)
148 );
149
150 labelList orderedToOld(sloan_order.size());
151 forAll(orderedToOld, c)
152 {
153 orderedToOld[c] = index_map[sloan_order[c]];
154 }
155
156 if (useReverse)
157 {
158 Foam::reverse(orderedToOld);
159 }
160
161 return orderedToOld;
163
164} // End anonymous namespace
165
166
167// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
168
170(
171 const polyMesh& mesh
172) const
173{
174 // Construct graph : faceOwner + connections across cyclics.
175
176 // Determine neighbour cell
177 labelList nbr(mesh.nBoundaryFaces(), -1);
178 for (const polyPatch& pp : mesh.boundaryMesh())
179 {
180 if (pp.coupled() && !isA<processorPolyPatch>(pp))
181 {
182 SubList<label>(nbr, pp.size(), pp.offset()) = pp.faceCells();
183 }
184 }
186
187
188 Graph G(mesh.nCells());
189
190 // Add internal faces
191 for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
192 {
193 add_edge(mesh.faceOwner()[facei], mesh.faceNeighbour()[facei], G);
194 }
195
196 // Add cyclics
197 for (const polyPatch& pp : mesh.boundaryMesh())
198 {
199 if
200 (
201 pp.coupled()
204 )
205 {
206 label bFacei = pp.offset();
207
208 for (const label ownCelli : pp.faceCells())
209 {
210 const label nbrCelli = nbr[bFacei];
211 ++bFacei;
212
213 if (ownCelli < nbrCelli)
214 {
215 add_edge(ownCelli, nbrCelli, G);
216 }
217 else
218 {
219 add_edge(nbrCelli, ownCelli, G);
220 }
221 }
223 }
224
225 return renumberImpl(G, reverse_);
226}
227
228
230(
231 const CompactListList<label>& cellCells
232) const
233{
234 Graph G(cellCells.size());
235
236 forAll(cellCells, celli)
237 {
238 const auto& neighbours = cellCells[celli];
239
240 for (const label nbr : neighbours)
241 {
242 if (celli < nbr)
243 {
244 add_edge(celli, nbr, G);
245 }
247 }
248
249 return renumberImpl(G, reverse_);
250}
251
252
254(
255 const labelListList& cellCells
256) const
257{
258 Graph G(cellCells.size());
259
260 forAll(cellCells, celli)
261 {
262 const auto& neighbours = cellCells[celli];
263
264 for (const label nbr : neighbours)
265 {
266 if (celli < nbr)
267 {
268 add_edge(celli, nbr, G);
269 }
270 }
271 }
272
273 return renumberImpl(G, reverse_);
274}
275
276
277// ************************************************************************* //
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
graph_traits< Graph >::vertices_size_type size_type
graph_traits< Graph >::vertex_descriptor Vertex
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A packed storage of objects of type <T> using an offset table for access.
label size() const noexcept
The primary size (the number of rows/sublists).
Sloan renumbering algorithm.
void reverse(bool on) noexcept
Toggle reverse on/off.
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.
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
label nInternalFaces() const noexcept
Number of internal faces.
label nCells() const noexcept
Number of mesh cells.
Abstract base class for renumbering.
renumberMethod()
Default construct.
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
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const dimensionedScalar c
Speed of light in a vacuum.
const dimensionedScalar G
Newtonian constant of gravitation.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition UListI.H:539
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
pointField vertices(const blockVertexList &bvl)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299