zoltanDecomp.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2021-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "zoltanDecomp.H"
27 #include "globalIndex.H"
28 #include "PstreamGlobals.H"
30 
31 #include "sigFpe.H"
32 #ifdef LINUX_GNUC
33  #ifndef __USE_GNU
34  #define __USE_GNU
35  #endif
36  #include <fenv.h>
37 #endif
38 
39 #include "zoltan.h"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
46 
48  (
51  distributor
52  );
53 }
54 
55 
56 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57 
58 static int get_number_of_vertices(void *data, int *ierr)
59 {
60  const Foam::pointField& points =
61  *static_cast<const Foam::pointField*>(data);
62 
63  *ierr = ZOLTAN_OK;
64 
65  return points.size();
66 }
67 
68 
69 static void get_vertex_list
70 (
71  void* data,
72  int nGID,
73  int nLID,
74  ZOLTAN_ID_PTR globalIDs,
75  ZOLTAN_ID_PTR localIDs,
76  int wgt_dim,
77  float* obj_wgts,
78  int* ierr
79 )
80 {
82  vertexData =
83  *static_cast
84  <
85  const Foam::Tuple2
86  <
87  const Foam::pointField&,
88  const Foam::scalarField&
89  >*
90  >(data);
91 
92  const Foam::pointField& points = vertexData.first();
93  const Foam::scalarField& weights = vertexData.second();
94 
95  if (wgt_dim != weights.size()/points.size())
96  {
97  *ierr = ZOLTAN_FATAL;
98  return;
99  }
100 
101  Foam::globalIndex globalMap(points.size());
102 
103  for (Foam::label i=0; i<points.size(); i++)
104  {
105  localIDs[i] = i;
106  globalIDs[i] = globalMap.toGlobal(i);
107 
108  for(int j=0; j<wgt_dim; j++)
109  {
110  obj_wgts[wgt_dim*i + j] = weights[wgt_dim*i + j];
111  }
112  }
113 
114  *ierr = ZOLTAN_OK;
115 }
116 
117 
118 static int get_mesh_dim(void* data, int* ierr)
119 {
120  return 3;
121 }
122 
123 
124 static void get_geom_list
125 (
126  void* data,
127  int nGID,
128  int nLID,
129  int nPoints,
130  ZOLTAN_ID_PTR globalIDs,
131  ZOLTAN_ID_PTR localIDs,
132  int nDim,
133  double* vertices,
134  int* ierr
135 )
136 {
137  const Foam::pointField& points =
138  *static_cast<const Foam::pointField*>(data);
139 
140  if
141  (
142  (nGID != 1)
143  || (nLID != 1)
144  || (nPoints != points.size())
145  || (nDim != 3)
146  )
147  {
148  *ierr = ZOLTAN_FATAL;
149  return;
150  }
151 
152  double* p = vertices;
153 
154  for (Foam::label celli=0; celli<nPoints; celli++)
155  {
156  const Foam::point& pt = points[celli];
157 
158  for (Foam::direction cmpt=0; cmpt<Foam::vector::nComponents; cmpt++)
159  {
160  *p++ = pt[cmpt];
161  }
162  }
163 
164  *ierr = ZOLTAN_OK;
165 }
166 
167 
169 (
170  void* data,
171  int nGID,
172  int nLID,
173  int nPoints,
174  ZOLTAN_ID_PTR globalIDs,
175  ZOLTAN_ID_PTR localIDs,
176  int* numEdges,
177  int* ierr
178 )
179 {
180  const Foam::CompactListList<Foam::label>& adjacency =
181  *static_cast<const Foam::CompactListList<Foam::label>*>(data);
182 
183  Foam::labelList sizes(adjacency.sizes());
184 
185  if ((nGID != 1) || (nLID != 1) || (nPoints != sizes.size()))
186  {
187  *ierr = ZOLTAN_FATAL;
188  return;
189  }
190 
191  for (Foam::label i=0; i<nPoints; i++)
192  {
193  numEdges[i] = sizes[i];
194  }
195 
196  *ierr = ZOLTAN_OK;
197 }
198 
199 
200 static void get_edge_list
201 (
202  void* data,
203  int nGID,
204  int nLID,
205  int nPoints,
206  ZOLTAN_ID_PTR globalIDs,
207  ZOLTAN_ID_PTR localIDs,
208  int* num_edges,
209  ZOLTAN_ID_PTR nborGID,
210  int* nborProc,
211  int wgt_dim,
212  float* ewgts,
213  int* ierr
214 )
215 {
216  const Foam::CompactListList<Foam::label>& adjacency =
217  *static_cast<const Foam::CompactListList<Foam::label>*>(data);
218 
219  const Foam::labelList& m(adjacency.m());
220 
221  Foam::globalIndex globalMap(nPoints);
222 
223  if
224  (
225  (nGID != 1)
226  || (nLID != 1)
227  )
228  {
229  *ierr = ZOLTAN_FATAL;
230  return;
231  }
232 
233  forAll(m, i)
234  {
235  nborGID[i] = m[i];
236  nborProc[i] = globalMap.whichProcID(m[i]);
237  }
238 
239  *ierr = ZOLTAN_OK;
240 }
241 
242 
244 (
245  const CompactListList<label>& adjacency,
246  const pointField& points,
247  const scalarField& pWeights,
248  List<label>& finalDecomp
249 ) const
250 {
251  if (Pstream::parRun() && !points.size())
252  {
253  Pout<< "No points on processor " << Pstream::myProcNo() << nl
254  << " Zoltan cannot redistribute without points"
255  "on all processors." << endl;
257  }
258 
259  stringList args(1);
260  args[0] = "zoltanDecomp";
261 
262  int argc = args.size();
263  char* argv[argc];
264  for (label i = 0; i < argc; i++)
265  {
266  argv[i] = strdup(args[i].c_str());
267  }
268 
269  float ver;
270  int rc = Zoltan_Initialize(argc, argv, &ver);
271 
272  if (rc != ZOLTAN_OK)
273  {
275  << "Failed initialising Zoltan" << exit(FatalError);
276  }
277 
278  struct Zoltan_Struct *zz = Zoltan_Create(PstreamGlobals::MPI_COMM_FOAM);
279 
280  const int nWeights = pWeights.size()/points.size();
281 
282  // Set internal parameters
283  Zoltan_Set_Param(zz, "return_lists", "export");
284  Zoltan_Set_Param(zz, "obj_weight_dim", name(nWeights).c_str());
285  Zoltan_Set_Param(zz, "edge_weight_dim", "0");
286 
287  // General default parameters
288  Zoltan_Set_Param(zz, "debug_level", "0");
289  Zoltan_Set_Param(zz, "imbalance_tol", "1.05");
290 
291  // Set default method parameters
292  Zoltan_Set_Param(zz, "lb_method", "graph");
293  Zoltan_Set_Param(zz, "lb_approach", "repartition");
294 
295  if (decompositionDict_.found("zoltanCoeffs"))
296  {
297  const dictionary& coeffsDict_ =
298  decompositionDict_.subDict("zoltanCoeffs");
299 
300  forAllConstIter(IDLList<entry>, coeffsDict_, iter)
301  {
302  if (!iter().isDict())
303  {
304  const word& key = iter().keyword();
305  const word value(iter().stream());
306 
307  if (debug)
308  {
309  Info<< typeName << " : setting parameter " << key
310  << " to " << value << endl;
311  }
312 
313  Zoltan_Set_Param(zz, key.c_str(), value.c_str());
314  }
315  }
316  }
317 
318  void* pointsPtr = &const_cast<pointField&>(points);
319  Zoltan_Set_Num_Obj_Fn(zz, get_number_of_vertices, pointsPtr);
320 
321  Tuple2<const pointField&, const scalarField&> vertexData(points, pWeights);
322  Zoltan_Set_Obj_List_Fn(zz, get_vertex_list, &vertexData);
323 
324  // Callbacks for geometry
325  Zoltan_Set_Num_Geom_Fn(zz, get_mesh_dim, pointsPtr);
326  Zoltan_Set_Geom_Multi_Fn(zz, get_geom_list, pointsPtr);
327 
328  // Callbacks for connectivity
329  void* adjacencyPtr = &const_cast<CompactListList<label>&>(adjacency);
330  Zoltan_Set_Num_Edges_Multi_Fn(zz, get_num_edges_list, adjacencyPtr);
331  Zoltan_Set_Edge_List_Multi_Fn(zz, get_edge_list, adjacencyPtr);
332 
333  int changed, num_imp, num_exp, *imp_procs, *exp_procs;
334  int *imp_to_part, *exp_to_part;
335  int num_gid_entries, num_lid_entries;
336  ZOLTAN_ID_PTR imp_global_ids, exp_global_ids;
337  ZOLTAN_ID_PTR imp_local_ids, exp_local_ids;
338 
339  // Switch off FPU error trapping to work around divide-by-zero bug
340  // in the Zoltan graph and hypergraph methods
341  #ifdef FE_NOMASK_ENV
342  int oldExcepts = fedisableexcept
343  (
344  FE_DIVBYZERO
345  | FE_INVALID
346  | FE_OVERFLOW
347  );
348  #endif
349 
350  // Perform load balancing
351  Zoltan_LB_Partition
352  (
353  zz,
354  &changed,
355  &num_gid_entries,
356  &num_lid_entries,
357  &num_imp,
358  &imp_global_ids,
359  &imp_local_ids,
360  &imp_procs,
361  &imp_to_part,
362  &num_exp,
363  &exp_global_ids,
364  &exp_local_ids,
365  &exp_procs,
366  &exp_to_part
367  );
368 
369  // Reinstate FPU error trapping
370  #ifdef FE_NOMASK_ENV
371  feenableexcept(oldExcepts);
372  #endif
373 
374  if (debug && changed)
375  {
376  Pout << "Zoltan number to move " << changed << " " << num_exp << endl;
377  }
378 
379  finalDecomp.setSize(points.size());
380  finalDecomp = Pstream::myProcNo();
381 
382  if (changed)
383  {
384  for(int i=0; i<num_exp; i++)
385  {
386  finalDecomp[exp_local_ids[i]] = exp_procs[i];
387  }
388  }
389 
390  // Free memory allocated for load-balancing results by Zoltan library
391 
392  Zoltan_LB_Free_Part
393  (
394  &imp_global_ids,
395  &imp_local_ids,
396  &imp_procs,
397  &imp_to_part
398  );
399 
400  Zoltan_LB_Free_Part
401  (
402  &exp_global_ids,
403  &exp_local_ids,
404  &exp_procs,
405  &exp_to_part
406  );
407 
408  Zoltan_Destroy(&zz);
409 
410  return 0;
411 }
412 
413 
414 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
415 
417 :
418  decompositionMethod(decompositionDict)
419 {}
420 
421 
422 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
423 
425 (
426  const polyMesh& mesh,
427  const pointField& cellCentres,
428  const scalarField& cellWeights
429 )
430 {
431  if (cellCentres.size() != mesh.nCells())
432  {
434  << "Can use this decomposition method only for the whole mesh"
435  << endl
436  << "and supply one coordinate (cellCentre) for every cell." << endl
437  << "The number of coordinates " << cellCentres.size() << endl
438  << "The number of cells in the mesh " << mesh.nCells()
439  << exit(FatalError);
440  }
441 
442 
443  // Convert mesh.cellCells() into compressed format
444  CompactListList<label> cellCells;
445  calcCellCells
446  (
447  mesh,
448  identityMap(mesh.nCells()),
449  mesh.nCells(),
450  true,
451  cellCells
452  );
453 
454  // Decompose using default weights
455  List<label> finalDecomp;
456  decompose
457  (
458  cellCells,
459  cellCentres,
460  cellWeights,
461  finalDecomp
462  );
463 
464  // Copy back to labelList
465  labelList decomp(cellCentres.size());
466  forAll(decomp, i)
467  {
468  decomp[i] = finalDecomp[i];
469  }
470  return decomp;
471 }
472 
473 
475 (
476  const polyMesh& mesh,
477  const labelList& agglom,
478  const pointField& agglomPoints,
479  const scalarField& pointWeights
480 )
481 {
482  if (agglom.size() != mesh.nCells())
483  {
485  << "Size of cell-to-coarse map " << agglom.size()
486  << " differs from number of cells in mesh " << mesh.nCells()
487  << exit(FatalError);
488  }
489 
490  // Convert agglom into compressed format
491  CompactListList<label> cellCells;
492  calcCellCells
493  (
494  mesh,
495  agglom,
496  agglomPoints.size(),
497  true,
498  cellCells
499  );
500 
501  // Decompose using weights
502  List<label> finalDecomp;
503  decompose
504  (
505  cellCells,
506  agglomPoints,
507  pointWeights,
508  finalDecomp
509  );
510 
511  // Rework back into decomposition for original mesh
512  labelList fineDistribution(agglom.size());
513 
514  forAll(fineDistribution, i)
515  {
516  fineDistribution[i] = finalDecomp[agglom[i]];
517  }
518 
519  return fineDistribution;
520 }
521 
522 
524 (
525  const labelListList& globalPointPoints,
526  const pointField& points,
527  const scalarField& pWeights
528 )
529 {
530  if (points.size() != globalPointPoints.size())
531  {
533  << "Inconsistent number of cells (" << globalPointPoints.size()
534  << ") and number of cell centres (" << points.size()
535  << ")." << exit(FatalError);
536  }
537 
538 
539  // Convert globalPointPoints into compressed format
540  CompactListList<label> pointPoints(globalPointPoints);
541 
542  // Decompose using weights
543  List<label> finalDecomp;
544  decompose
545  (
546  pointPoints,
547  points,
548  pWeights,
549  finalDecomp
550  );
551 
552  return finalDecomp;
553 }
554 
555 
556 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
A packed storage unstructured matrix of objects of type <T> using an offset table for access.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:63
const Type2 & second() const
Return second.
Definition: Tuple2.H:128
const Type1 & first() const
Return first.
Definition: Tuple2.H:116
const UList< T > & m() const
Return the packed matrix of data.
label size() const
Return the primary size, i.e. the number of rows.
labelList sizes() const
Return sizes (to be used e.g. for construction)
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
label size() const
Return the number of arguments.
Definition: argListI.H:90
Abstract base class for decomposition.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
label whichProcID(const label i) const
Which processor does global come from? Binary search.
Definition: globalIndexI.H:123
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
label nCells() const
Zoltan redistribution in parallel.
Definition: zoltanDecomp.H:95
virtual labelList decompose(const pointField &points, const scalarField &pointWeights)
Inherit decompose from decompositionMethod.
zoltanDecomp(const dictionary &decompositionDict)
Construct given the decomposition dictionary and mesh.
Definition: zoltanDecomp.C:416
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
volScalarField scalarField(fieldObject, mesh)
const pointField & points
label nPoints
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
pointField vertices(const blockVertexList &bvl)
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:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
List< string > stringList
A List of strings.
Definition: stringList.H:50
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
uint8_t direction
Definition: direction.H:45
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::argList args(argc, argv)
volScalarField & p
static int get_number_of_vertices(void *data, int *ierr)
Definition: zoltanDecomp.C:58
static void get_num_edges_list(void *data, int nGID, int nLID, int nPoints, ZOLTAN_ID_PTR globalIDs, ZOLTAN_ID_PTR localIDs, int *numEdges, int *ierr)
Definition: zoltanDecomp.C:169
static void get_geom_list(void *data, int nGID, int nLID, int nPoints, ZOLTAN_ID_PTR globalIDs, ZOLTAN_ID_PTR localIDs, int nDim, double *vertices, int *ierr)
Definition: zoltanDecomp.C:125
static int get_mesh_dim(void *data, int *ierr)
Definition: zoltanDecomp.C:118
static void get_vertex_list(void *data, int nGID, int nLID, ZOLTAN_ID_PTR globalIDs, ZOLTAN_ID_PTR localIDs, int wgt_dim, float *obj_wgts, int *ierr)
Definition: zoltanDecomp.C:70
static void get_edge_list(void *data, int nGID, int nLID, int nPoints, ZOLTAN_ID_PTR globalIDs, ZOLTAN_ID_PTR localIDs, int *num_edges, ZOLTAN_ID_PTR nborGID, int *nborProc, int wgt_dim, float *ewgts, int *ierr)
Definition: zoltanDecomp.C:201