zoltanRenumber.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) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::zoltanRenumber
28
29Description
30 Renumber using Zoltan
31
32 Zoltan install:
33 - in your ~/.bashrc:
34 export ZOLTAN_ARCH_DIR=\
35 $WM_THIRD_PARTY_DIR/platforms/linux64Gcc/Zoltan_XXX
36 - unpack into $WM_THIRD_PARTY_DIR
37 - cd Zoltan_XXX
38 - mkdir build
39 - cd build
40 - export CCFLAGS="-fPIC"
41 - export CXXFLAGS="-fPIC"
42 - export CFLAGS="-fPIC"
43 - export LDFLAGS="-shared"
44 - ../configure \
45 --prefix=$ZOLTAN_ARCH_DIR \
46 --with-ccflags=-fPIC --with-cxxflags=-fPIC --with-ldflags=-shared
47
48SourceFiles
49 zoltanRenumber.C
50
51\*---------------------------------------------------------------------------*/
52
53#include "zoltanRenumber.H"
55#include "IFstream.H"
56#include "labelIOList.H"
57#include "polyMesh.H"
58#include "globalMeshData.H"
59#include "globalIndex.H"
60#include "uint.H"
61
62#pragma GCC diagnostic ignored "-Wold-style-cast"
63#include "zoltan.h"
64#include <mpi.h>
65
66// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67
68namespace Foam
73 (
77 );
78}
79
81static int get_number_of_vertices(void *data, int *ierr)
82{
83 const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
84 *ierr = ZOLTAN_OK;
85 return mesh.nCells();
86}
87
89static void get_vertex_list(void *data, int sizeGID, int sizeLID,
90 ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID,
91 int wgt_dim, float *obj_wgts, int *ierr)
92{
93 const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
94
95 *ierr = ZOLTAN_OK;
96
97 /* In this example, return the IDs of our vertices, but no weights.
98 * Zoltan will assume equally weighted vertices.
99 */
100
101 wgt_dim = 0;
102 obj_wgts = nullptr;
103
104 for (Foam::label i=0; i<mesh.nCells(); i++)
105 {
106 globalID[i] = i; // should be global
107 localID[i] = i;
108 }
109}
110
112static void get_num_edges_list(void *data, int sizeGID, int sizeLID,
113 int num_obj,
114 ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID,
115 int *numEdges, int *ierr)
116{
117 const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
118
119 if ((sizeGID != 1) || (sizeLID != 1) || (num_obj != mesh.nCells()))
120 {
121 *ierr = ZOLTAN_FATAL;
122 return;
123 }
124
125 for (Foam::label i=0; i < num_obj ;i++)
126 {
127 Foam::label celli = localID[i];
128 const Foam::cell& cFaces = mesh.cells()[celli];
129 forAll(cFaces, cFacei)
130 {
131 Foam::label n = 0;
132 if (mesh.isInternalFace(cFaces[cFacei]))
133 {
134 n++;
135 }
136 numEdges[i] = n;
137 }
138 }
139
140 *ierr = ZOLTAN_OK;
141}
142
144static void get_edge_list(void *data, int sizeGID, int sizeLID,
145 int num_obj, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID,
146 int *num_edges,
147 ZOLTAN_ID_PTR nborGID, int *nborProc,
148 int wgt_dim, float *ewgts, int *ierr)
149{
150 const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
151
152 if
153 (
154 (sizeGID != 1)
155 || (sizeLID != 1)
156 || (num_obj != mesh.nCells())
157 || (wgt_dim != 1)
158 )
159 {
160 *ierr = ZOLTAN_FATAL;
161 return;
162 }
163
164 ZOLTAN_ID_TYPE* nextNbor = nborGID;
165 int* nextProc = nborProc;
166 float* nextWgt = ewgts;
167
168 for (Foam::label i=0; i < num_obj; i++)
169 {
170 Foam::label celli = localID[i];
171
172 const Foam::cell& cFaces = mesh.cells()[celli];
173 forAll(cFaces, cFacei)
174 {
175 Foam::label n = 0;
176
177 Foam::label facei = cFaces[cFacei];
178 if (mesh.isInternalFace(facei))
179 {
180 Foam::label nbr = mesh.faceOwner()[facei];
181 if (nbr == celli)
182 {
183 nbr = mesh.faceNeighbour()[facei];
184 }
185
186 // Note: global index
187 *nextNbor++ = nbr;
188 *nextProc++ = 0;
189 *nextWgt++ = 1.0;
190
191 n++;
192 }
193 if (n != num_edges[i])
194 {
195 *ierr = ZOLTAN_FATAL;
196 return;
197 }
198 }
199 }
200 *ierr = ZOLTAN_OK;
201}
202
204static int get_mesh_dim(void *data, int *ierr)
205{
206 const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
207
208 return mesh.nSolutionD();
209}
210
212static void get_geom_list
213(
214 void *data,
215 int num_gid_entries,
216 int num_lid_entries,
217 int num_obj,
218 ZOLTAN_ID_PTR global_ids,
219 ZOLTAN_ID_PTR local_ids,
220 int num_dim,
221 double *geom_vec,
222 int *ierr
223)
224{
225 const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
226
227 if
228 (
229 (num_gid_entries != 1)
230 || (num_lid_entries != 1)
231 || (num_obj != mesh.nCells())
232 || (num_dim != mesh.nSolutionD())
233 )
234 {
235 *ierr = ZOLTAN_FATAL;
236 return;
237 }
238
239 double* p = geom_vec;
240
241
242 const Foam::Vector<Foam::label>& sol = mesh.solutionD();
243
244 const Foam::pointField& cc = mesh.cellCentres();
245
246 for (Foam::label celli = 0; celli < num_obj; celli++)
247 {
248 const Foam::point& pt = cc[celli];
249
250 for (Foam::direction cmpt = 0; cmpt < Foam::vector::nComponents; cmpt++)
251 {
252 if (sol[cmpt] == 1)
253 {
254 *p++ = pt[cmpt];
255 }
256 }
257 }
258 *ierr = ZOLTAN_OK;
259}
260
261
262// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
265:
267 coeffsDict_(dict.optionalSubDict(typeName+"Coeffs"))
268{}
269
270
271// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
274(
275 const polyMesh& pMesh,
276 const pointField& points
277) const
278{
279 stringList args(1);
280 args[0] = "zoltanRenumber";
281
282 int argc = args.size();
283 char* argv[argc];
284 for (label i = 0; i < argc; i++)
285 {
286 argv[i] = strdup(args[i].c_str());
287 }
288
289 float ver;
290 int rc = Zoltan_Initialize(argc, argv, &ver);
291
292 Foam::Pout<< "Initialised to " << ver << Foam::endl;
293
294 if (rc != ZOLTAN_OK)
295 {
297 << "Failed initialising Zoltan" << exit(FatalError);
298 }
299
300 struct Zoltan_Struct *zz = Zoltan_Create(MPI_COMM_WORLD);
301
302 polyMesh& mesh = const_cast<polyMesh&>(pMesh);
303
304
305 for (const entry& dEntry : coeffsDict_)
306 {
307 if (!dEntry.isDict())
308 {
309 const word& key = dEntry.keyword();
310 const word value(dEntry.get<word>());
311
312 Info<< typeName << " : setting parameter " << key
313 << " to " << value << endl;
314
315 Zoltan_Set_Param(zz, key.c_str(), value.c_str());
316 }
317 }
318
319
320 Zoltan_Set_Num_Obj_Fn(zz, get_number_of_vertices, &mesh);
321 Zoltan_Set_Obj_List_Fn(zz, get_vertex_list, &mesh);
322
323 // Callbacks for geometry
324 Zoltan_Set_Num_Geom_Fn(zz, get_mesh_dim, &mesh);
325 Zoltan_Set_Geom_Multi_Fn(zz, get_geom_list, &mesh);
326
327 // Callbacks for connectivity
328 Zoltan_Set_Num_Edges_Multi_Fn(zz, get_num_edges_list, &mesh);
329 Zoltan_Set_Edge_List_Multi_Fn(zz, get_edge_list, &mesh);
330
331
332
333 //Note: !global indices
334 List<ZOLTAN_ID_TYPE> wantedCells(mesh.nCells());
335
336 globalIndex globalCells(mesh.nCells());
337 forAll(wantedCells, i)
338 {
339 //wantedCells[i] = i;
340 wantedCells[i] = globalCells.toGlobal(i);
341 }
342
343 List<ZOLTAN_ID_TYPE> oldToNew(mesh.nCells());
344
345 int err = Zoltan_Order
346 (
347 zz,
348 1, //int num_gid_entries,
349 mesh.globalData().nTotalCells(), //int num_obj,
350 wantedCells.begin(),
351 oldToNew.begin()
352 );
353
354 if (err != ZOLTAN_OK)
355 {
357 << "Failed Zoltan_Order" << exit(FatalError);
358 }
359
360
361 for (label i = 0; i < argc; i++)
362 {
363 free(argv[i]);
364 }
365
366
367 labelList order(oldToNew.size());
368 forAll(order, i)
369 {
370 order[i] = oldToNew[i];
371 }
372 return order;
373}
374
375
376// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:329
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:65
label size() const noexcept
The number of arguments.
Definition: argListI.H:146
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:70
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
label toGlobal(const label i) const
From local to global index.
Definition: globalIndexI.H:260
label nTotalCells() const noexcept
Return total number of cells in decomposed mesh.
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1310
label nCells() const noexcept
Number of mesh cells.
Abstract base class for renumbering.
A class for handling words, derived from Foam::string.
Definition: word.H:68
Renumber using Zoltan.
virtual labelList renumber(const pointField &) const
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
uint8_t direction
Definition: direction.H:56
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dictionary dict
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
System unsigned integer.
static void get_geom_list(void *data, int num_gid_entries, int num_lid_entries, int num_obj, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int num_dim, double *geom_vec, int *ierr)
static int get_number_of_vertices(void *data, int *ierr)
static int get_mesh_dim(void *data, int *ierr)
static void get_edge_list(void *data, int sizeGID, int sizeLID, int num_obj, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int *num_edges, ZOLTAN_ID_PTR nborGID, int *nborProc, int wgt_dim, float *ewgts, int *ierr)
static void get_num_edges_list(void *data, int sizeGID, int sizeLID, int num_obj, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int *numEdges, int *ierr)
static void get_vertex_list(void *data, int sizeGID, int sizeLID, ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID, int wgt_dim, float *obj_wgts, int *ierr)