Loading...
Searching...
No Matches
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 Copyright (C) 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
27Notes
28 Renumber using Zoltan
29
30 Zoltan install:
31 - in your ~/.bashrc:
32 export ZOLTAN_ARCH_DIR=\
33 $WM_THIRD_PARTY_DIR/platforms/linux64Gcc/Zoltan_XXX
34 - unpack into $WM_THIRD_PARTY_DIR
35 - cd Zoltan_XXX
36 - mkdir build
37 - cd build
38 - export CCFLAGS="-fPIC"
39 - export CXXFLAGS="-fPIC"
40 - export CFLAGS="-fPIC"
41 - export LDFLAGS="-shared"
42 - ../configure \
43 --prefix=$ZOLTAN_ARCH_DIR \
44 --with-ccflags=-fPIC --with-cxxflags=-fPIC --with-ldflags=-shared
45
46\*---------------------------------------------------------------------------*/
47
48#include "zoltanRenumber.H"
50#include "IFstream.H"
51#include "labelIOList.H"
52#include "polyMesh.H"
53#include "globalMeshData.H"
54#include "globalIndex.H"
55#include "CStringList.H"
56#include "uint.H"
57#include <algorithm>
58#include <numeric>
59
60// Include MPI without any C++ bindings
61#ifndef MPICH_SKIP_MPICXX
62#define MPICH_SKIP_MPICXX
63#endif
64#ifndef OMPI_SKIP_MPICXX
65#define OMPI_SKIP_MPICXX
66#endif
67
68#pragma GCC diagnostic ignored "-Wold-style-cast"
69#include "zoltan.h"
70#include <mpi.h>
72// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
73
74namespace Foam
75{
77
79 (
83 );
85
86
87// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
88
89// [ZOLTAN_NUM_OBJ_FN]
90// The number of graph vertices (locally)
91static int get_number_of_vertices(void *data, int *ierr)
92{
93 *ierr = ZOLTAN_OK;
94 const auto& mesh = *static_cast<const Foam::polyMesh*>(data);
96 return mesh.nCells();
97}
98
99
100// [ZOLTAN_OBJ_LIST_FN]
101// The ids for the graph vertices
102static void get_vertex_list
103(
104 void *data,
105 int sizeGID, // (currently unused)
106 int sizeLID, // (currently unused)
107 ZOLTAN_ID_PTR global_ids,
108 ZOLTAN_ID_PTR local_ids,
109 int wgt_dim, // (currently unused)
110 float *obj_wgts, // (currently unused)
111 int *ierr
112)
113{
114 *ierr = ZOLTAN_OK;
115 const auto& mesh = *static_cast<const Foam::polyMesh*>(data);
116
117 const auto nCells = mesh.nCells();
118
119 // Offset for globally unique IDs
120 const auto myProci = Foam::UPstream::myProcNo();
121 const auto myProcOffset =
122 mesh.globalData().globalMeshCellAddr().localStart(myProci);
123
124 // Global indices
125 std::iota(global_ids, (global_ids + nCells), ZOLTAN_ID_TYPE(myProcOffset));
126
127 // Local indices
128 std::iota(local_ids, (local_ids + nCells), ZOLTAN_ID_TYPE(0));
129
130 // No weights?
131 // Zoltan will assume equally weighted vertices.
132
133 // if (obj_wgts && wgt_dim > 0)
134 // {
135 // std::fill_n(obj_wgts, wgt_dim, float(1));
136 // }
137}
138
139
140// [ZOLTAN_NUM_EDGES_MULTI_FN]
141// query function returns the number of edges in the communication
142// graph of the application for each object in a list of objects.
143// That is, for each object in the global_ids/local_ids arrays, the
144// number of objects with which the given object must share
145// information is returned.
146
147static void get_num_edges_list
148(
149 void *data,
150 int sizeGID,
151 int sizeLID,
152 int num_obj, // Number of graph vertices (mesh cells)
153 ZOLTAN_ID_PTR global_ids, // (currently unused)
154 ZOLTAN_ID_PTR local_ids, // rank-local vertex id (mesh cell id)
155 int *numEdges,
156 int *ierr
157)
158{
159 *ierr = ZOLTAN_OK;
160 const auto& mesh = *static_cast<const Foam::polyMesh*>(data);
161
162 if ((sizeGID != 1) || (sizeLID != 1) || (num_obj != mesh.nCells()))
163 {
164 *ierr = ZOLTAN_FATAL;
165 return;
166 }
167
168 const Foam::label nCells = num_obj;
169
170 for (Foam::label i=0; i < nCells; ++i)
171 {
172 const Foam::label celli = local_ids[i];
173 int numNbr = 0;
174
175 for (const auto facei : mesh.cells()[celli])
176 {
177 if (mesh.isInternalFace(facei))
178 {
179 ++numNbr;
180 }
181 // TBD: check coupled etc
182 }
184 numEdges[i] = numNbr;
185 }
186}
187
188
189// [ZOLTAN_EDGE_LIST_MULTI_FN]
190static void get_edge_list
191(
192 void *data,
193 int sizeGID,
194 int sizeLID,
195 int num_obj,
196 ZOLTAN_ID_PTR global_ids, // (currently unused)
197 ZOLTAN_ID_PTR local_ids, // rank-local vertex id (mesh cell id)
198 int *num_edges,
199 ZOLTAN_ID_PTR nborGID,
200 int *nborProc,
201 int wgt_dim,
202 float *ewgts,
203 int *ierr
204)
205{
206 *ierr = ZOLTAN_OK;
207 const auto& mesh = *static_cast<const Foam::polyMesh*>(data);
208
209 if
210 (
211 (sizeGID != 1)
212 || (sizeLID != 1)
213 || (num_obj != mesh.nCells())
214 || (wgt_dim != 1)
215 )
216 {
217 *ierr = ZOLTAN_FATAL;
218 return;
219 }
220
221 // Offset for globally unique IDs
222 const auto myProci = Foam::UPstream::myProcNo();
223 const auto myProcOffset =
224 mesh.globalData().globalMeshCellAddr().localStart(myProci);
225
226 auto* nextNbor = nborGID;
227 auto* nextProc = nborProc;
228 auto* nextWgt = ewgts;
229
230 const Foam::label nCells = num_obj;
231
232 for (Foam::label i=0; i < nCells; ++i)
233 {
234 const Foam::label celli = local_ids[i];
235 int numNbr = 0;
236
237 for (const auto facei : mesh.cells()[celli])
238 {
239 if (mesh.isInternalFace(facei))
240 {
241 Foam::label nbr = mesh.faceOwner()[facei];
242 if (nbr == celli)
243 {
244 nbr = mesh.faceNeighbour()[facei];
245 }
246
247 *nextNbor++ = (nbr + myProcOffset); // global id
248 *nextProc++ = myProci; // rank-local connection
249 *nextWgt++ = 1.0;
250
251 ++numNbr;
252 }
253 }
254
255 // Sanity check
256 if (numNbr != num_edges[i])
257 {
258 *ierr = ZOLTAN_FATAL;
259 return;
261 }
262}
263
264
265// [ZOLTAN_NUM_GEOM_FN]
266// The dimensionality of the mesh geometry (2D,3D, etc)
267static int get_mesh_dim(void *data, int *ierr)
268{
269 *ierr = ZOLTAN_OK;
270 const auto& mesh = *static_cast<const Foam::polyMesh*>(data);
272 return mesh.nSolutionD();
273}
274
275
276// [ZOLTAN_GEOM_MULTI_FN]
277// The geometric location of the graph vertices (mesh cellCentres)
278static void get_geom_list
279(
280 void *data,
281 int num_gid_entries,
282 int num_lid_entries,
283 int num_obj,
284 ZOLTAN_ID_PTR global_ids,
285 ZOLTAN_ID_PTR local_ids,
286 int num_dim, // dimensionality (2|3)
287 double *geom_vec, // [out] cellCentres
288 int *ierr
289)
290{
291 *ierr = ZOLTAN_OK;
292 const Foam::polyMesh& mesh = *static_cast<const Foam::polyMesh*>(data);
293
294 if
295 (
296 (num_gid_entries != 1)
297 || (num_lid_entries != 1)
298 || (num_obj != mesh.nCells())
299 || (num_dim != mesh.nSolutionD())
300 )
301 {
302 *ierr = ZOLTAN_FATAL;
303 return;
304 }
305
306 const Foam::pointField& cc = mesh.cellCentres();
307
308 if (num_dim == 3)
309 {
310 // Fast path for 3D - assumes that local_ids are still ordered
311 std::copy_n
312 (
313 reinterpret_cast<const Foam::point::cmptType*>(cc.cdata()),
314 (3*num_obj),
315 geom_vec
316 );
317 }
318 else
319 {
320 // [out] cellCentres
321 double* p = geom_vec;
322
323 const Foam::Vector<Foam::label>& sol = mesh.solutionD();
324
325 const Foam::label nCells = num_obj;
326
327 // Assumes that local_ids are still ordered
328 for (Foam::label celli = 0; celli < nCells; ++celli)
329 {
330 const Foam::point& pt = cc[celli];
331
332 for
333 (
334 Foam::direction cmpt = 0;
336 ++cmpt
337 )
338 {
339 if (sol[cmpt] == 1)
340 {
341 *p++ = pt[cmpt];
342 }
343 }
345 }
346}
347
348
349// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
350
351Foam::zoltanRenumber::zoltanRenumber(const dictionary& dict)
352:
354 coeffsDict_(dict.optionalSubDict(typeName + "Coeffs"))
355{}
356
357
358// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
359
361(
362 const polyMesh& pMesh
363) const
364{
365 // Zoltan_Initialize will trigger MPI_Init() if not already done.
366 // - use UPstream::initNull() so that OpenFOAM also knows about MPI
367
369
370 CStringList args({"zoltan-renumber"});
371
372 float ver;
373 int rc = Zoltan_Initialize(args.size(), args.strings(), &ver);
374
375 if (rc != ZOLTAN_OK)
376 {
378 << "Failed initialising Zoltan" << exit(FatalError);
379 }
380
381 struct Zoltan_Struct *zz = Zoltan_Create(MPI_COMM_WORLD);
382
383 // const bool verbose = coeffsDict_.getOrDefault(verbose, false);
384 const bool verbose = true;
385
386 // Default order method
387 Zoltan_Set_Param(zz, "ORDER_METHOD", "LOCAL_HSFC");
388
389 if (false)
390 {
391 Info<< typeName << " : default ORDER_METHOD = LOCAL_HSFC" << nl
392 << typeName << " : default ORDER_TYPE = LOCAL" << nl;
393 }
394
395 for (const entry& dEntry : coeffsDict_)
396 {
397 const word& key = dEntry.keyword();
398
399 // Internal keywords
400 if
401 (
402 key == "method"
403 || key == "verbose"
404 )
405 {
406 continue;
407 }
408
409 if (!dEntry.isDict())
410 {
411 const word value(dEntry.get<word>());
412
413 if (verbose)
414 {
415 Info<< typeName
416 << " : setting parameter "
417 << key << " = " << value << nl;
418 }
419
420 Zoltan_Set_Param(zz, key.c_str(), value.c_str());
421 }
422 }
423
424 // Always use rank LOCAL ordering
425 Zoltan_Set_Param(zz, "ORDER_TYPE", "LOCAL");
426
427
428 // Set callbacks
429 polyMesh& mesh = const_cast<polyMesh&>(pMesh);
430
431 // Callbacks for graph vertex IDs
432 Zoltan_Set_Num_Obj_Fn(zz, get_number_of_vertices, &mesh);
433 Zoltan_Set_Obj_List_Fn(zz, get_vertex_list, &mesh);
434
435 // Callbacks for geometry
436 Zoltan_Set_Num_Geom_Fn(zz, get_mesh_dim, &mesh);
437 Zoltan_Set_Geom_Multi_Fn(zz, get_geom_list, &mesh);
438
439 // Callbacks for connectivity
440 Zoltan_Set_Num_Edges_Multi_Fn(zz, get_num_edges_list, &mesh);
441 Zoltan_Set_Edge_List_Multi_Fn(zz, get_edge_list, &mesh);
442
443
444 const auto nCells = mesh.nCells();
445
446 List<ZOLTAN_ID_TYPE> globalIds(nCells);
447 List<ZOLTAN_ID_TYPE> oldToNew(nCells);
448
449 // Offset for globally unique IDs
450 const label myProci = UPstream::myProcNo();
451 const label myProcOffset =
452 mesh.globalData().globalMeshCellAddr().localStart(myProci);
453
454
455 // Global indices
456 std::iota(globalIds.begin(), globalIds.end(), ZOLTAN_ID_TYPE(myProcOffset));
457
458 int err = Zoltan_Order
459 (
460 zz,
461 1, //int num_gid_entries,
462 nCells, //int num_obj,
463 globalIds.begin(),
464 oldToNew.begin()
465 );
466
467 if (err != ZOLTAN_OK)
468 {
470 << "Failed Zoltan_Order" << exit(FatalError);
471 }
472
473 Zoltan_Destroy(&zz);
474
475
476 // From global number to rank-local numbering
477 globalIds.clear();
478 labelList order(oldToNew.size());
479
480 // From global to local (without checks)
481 std::transform
482 (
483 oldToNew.begin(),
484 oldToNew.end(),
485 order.begin(),
486 [=](const ZOLTAN_ID_TYPE id) -> label { return (id - myProcOffset); }
487 );
488
489 return order;
490}
491
492
493// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
An adapter for copying a list of C++ strings into a list of C-style strings for passing to C code tha...
Definition CStringList.H:68
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition UListI.H:410
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition UListI.H:454
const T * cdata() const noexcept
Return pointer to the underlying array serving as data storage.
Definition UListI.H:267
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool initNull()
Special purpose initialisation function.
Definition UPstream.C:30
static constexpr direction nComponents
Number of components in this vector space.
Cmpt cmptType
Component type.
Definition VectorSpace.H:91
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Abstract base class for renumbering.
renumberMethod()
Default construct.
A class for handling words, derived from Foam::string.
Definition word.H:66
virtual labelList renumber(const polyMesh &mesh) const
Return the cell visit order (from ordered back to original cell id) uses the mesh for connectivity an...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
Foam::argList args(argc, argv)
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 void get_vertex_list(void *data, int sizeGID, int sizeLID, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int wgt_dim, float *obj_wgts, int *ierr)
static int get_number_of_vertices(void *data, int *ierr)
static int get_mesh_dim(void *data, int *ierr)
static void get_num_edges_list(void *data, int sizeGID, int sizeLID, int num_obj, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int *numEdges, int *ierr)
static void get_edge_list(void *data, int sizeGID, int sizeLID, int num_obj, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids, int *num_edges, ZOLTAN_ID_PTR nborGID, int *nborProc, int wgt_dim, float *ewgts, int *ierr)