Arcane  4.1.12.0
Developer documentation
Loading...
Searching...
No Matches
ZoltanMeshPartitioner.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2026 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* ZoltanMeshPartitioner.cc (C) 2000-2025 */
9/* */
10/* Mesh partitioner using the Zoltan library. */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
14#include "arcane/utils/PlatformUtils.h"
15#include "arcane/utils/Convert.h"
16#include "arcane/utils/CheckedConvert.h"
17#include "arcane/utils/Array.h"
18#include "arcane/utils/StringBuilder.h"
19#include "arcane/utils/Iostream.h"
20#include "arcane/utils/ScopedPtr.h"
21#include "arcane/utils/NotImplementedException.h"
22#include "arcane/utils/ArgumentException.h"
23
24#include "arcane/core/ISubDomain.h"
25#include "arcane/core/IParallelMng.h"
27#include "arcane/core/IMeshSubMeshTransition.h"
28#include "arcane/core/ItemGroup.h"
29#include "arcane/core/Service.h"
30#include "arcane/core/Timer.h"
31#include "arcane/core/FactoryService.h"
32#include "arcane/core/ItemPrinter.h"
33#include "arcane/core/IItemFamily.h"
34#include "arcane/core/MeshVariable.h"
35#include "arcane/core/VariableBuildInfo.h"
36#include "arcane/core/CommonVariables.h"
37#include "arcane/utils/HashTableMap.h"
38
39#include "arcane_internal_config.h"
40
41// In case we use mpich2 or openmpi
42#define MPICH_SKIP_MPICXX
43#define OMPI_SKIP_MPICXX
44#include <zoltan.h>
45
46#include "arcane/std/MeshPartitionerBase.h"
47#include "arcane/std/ZoltanMeshPartitioner_axl.h"
48
49#include <set>
50
51#ifdef WIN32
52#include <Windows.h>
53#define sleep Sleep
54#endif
55
56#define ARCANE_DEBUG_ZOLTAN
57
58// GG: NOTE: The current implementation only supports uniqueId() on 32 bits.
59
60/*---------------------------------------------------------------------------*/
61/*---------------------------------------------------------------------------*/
62
63namespace Arcane
64{
65
66/*---------------------------------------------------------------------------*/
67/*---------------------------------------------------------------------------*/
68
69enum ZoltanModel
70{
71 ZOLTAN_HG = (1 << 0),
72 ZOLTAN_GRAPH = (1 << 1),
73 ZOLTAN_GEOM = (1 << 2)
74};
75
76/*---------------------------------------------------------------------------*/
77/*---------------------------------------------------------------------------*/
78
83class ZoltanInfo
84: public TraceAccessor
85{
86 public:
87
88 ZoltanInfo(MeshPartitionerBase* basePartitionner, ostream* ofile, int model = 1, const Real edgeWeightMultiplier = 1.)
89 : TraceAccessor(basePartitionner->mesh()->traceMng())
90 , m_mesh_partitionner_base(basePartitionner)
91 , m_nbEdges(0)
92 , m_nbPins(0)
93 , m_ofile(ofile)
94 , m_model(model)
95 , m_edgeGIDStart(0)
96 , m_edge_weight_multiplier(edgeWeightMultiplier)
97 {
98 m_own_cells = m_mesh_partitionner_base->mesh()->ownCells();
99 build();
100 }
101 MeshPartitionerBase* m_mesh_partitionner_base;
102
103 private:
104
105 CellGroup m_own_cells;
106 int m_nbEdges;
107 int m_nbPins;
108 ostream* m_ofile;
109 int m_model;
110 int m_edgeGIDStart;
111 Real m_edge_weight_multiplier;
112 std::set<std::pair<Int64, Int64>> m_weight_set;
113
114 public:
115
116 void build()
117 {
118 info() << "ZoltanInfo::build()";
119 if (m_model & ZOLTAN_GEOM) // No topological informations to build
120 return;
121
122 Integer nbOwnCells = m_own_cells.size();
123 Integer nbOwnEdges = 0;
124
125 if (m_model & ZOLTAN_HG) {
126 nbOwnEdges += nbOwnCells;
127 }
128 if (m_model & ZOLTAN_GRAPH) {
129 // Neighboors, each and then neighboorhood HE
130 nbOwnEdges += nbOwnCells * 6;
131 }
132
133 // Compute m_edgeGIDStart = Sum_k<i nbOwnEdges
134 IParallelMng* pm = m_mesh_partitionner_base->mesh()->parallelMng();
135 Int32UniqueArray scanouille(pm->commSize());
136 scanouille.fill(0);
137 scanouille[pm->commRank()] = nbOwnEdges;
138 pm->scan(Parallel::ReduceSum, scanouille);
139 m_edgeGIDStart = 0;
140 for (int i = 0; i < pm->commRank(); ++i) {
141 m_edgeGIDStart += scanouille[i];
142 }
143
144 m_nbEdges = 0;
145 m_nbPins = 0;
146 ENUMERATE_CELL (iCell, m_own_cells) {
147 if (!m_mesh_partitionner_base->cellUsedWithConstraints(*iCell))
148 continue;
149
150 int nbNgb = m_mesh_partitionner_base->nbNeighbourCellsWithConstraints(*iCell);
151 if (m_model & ZOLTAN_HG) {
152 m_nbEdges++;
153 m_nbPins += nbNgb + 1;
154 }
155 if (m_model & ZOLTAN_GRAPH) {
156 m_nbEdges += nbNgb;
157 m_nbPins += (nbNgb * 2);
158 }
159 }
160
161 info() << "nbEdges=" << m_nbEdges << " ; nbPins=" << m_nbPins;
162
163 if (m_mesh_partitionner_base->haveWeakConstraints()) {
164 MeshVariableScalarRefT<Face, Integer> weak_constraint(VariableBuildInfo(m_mesh_partitionner_base->mesh(), "EdgeWeight"));
165 ENUMERATE_FACE (iface, m_mesh_partitionner_base->mesh()->ownFaces()) {
166 const Face& face = *iface;
167 const Cell& bCell = iface->backCell();
168 const Cell& fCell = iface->frontCell();
169 if (bCell.null() || fCell.null())
170 continue;
171 if (weak_constraint[face] == 2) {
172 m_weight_set.insert(std::pair<Int64, Int64>(face.backCell().uniqueId(), face.frontCell().uniqueId()));
173 m_weight_set.insert(std::pair<Int64, Int64>(face.frontCell().uniqueId(), face.backCell().uniqueId()));
174 }
175 }
176 }
177 } // end build()
178
179 public:
180
181 static int getHgNumVertices(void* data, int* ierr)
182 {
183 ARCANE_UNUSED(ierr);
184
185 ZoltanInfo* zi = (ZoltanInfo*)data;
186
187 /*
188 * Supply this query function to Zoltan with Zoltan_Set_Num_Obj_Fn().
189 * It returns the number of vertices that this process owns.
190 *
191 * The parallel hypergraph method requires that vertex global IDs and
192 * weights are returned by the application with query functions. The
193 * global IDs should be unique, and no two processes should
194 * return the same vertex.
195 */
196 return zi->m_mesh_partitionner_base->nbOwnCellsWithConstraints();
197 }
198
199 static void getHgVerticesAndWeights(void* data, int num_gid_entries,
200 int num_lid_entries, ZOLTAN_ID_PTR gids, ZOLTAN_ID_PTR lids,
201 int wgt_dim, float* obj_weights, int* ierr)
202 {
203 ARCANE_UNUSED(num_gid_entries);
204 ARCANE_UNUSED(num_lid_entries);
205
206 ZoltanInfo* zi = (ZoltanInfo*)data;
207 /*
208 * Supply this query function to Zoltan with Zoltan_Set_Obj_List_Fn().
209 *
210 * It supplies vertex global IDs, local IDs,
211 * and weights to the Zoltan library. The application has previously
212 * indicated the number of weights per vertex by setting the
213 * OBJ_WEIGHT_DIM parameter.
214 */
215 //int nb_cell = zi->m_own_cells.size();
216
217 SharedArray<float> cells_weights;
218 if (wgt_dim != 0) { // No weight, does not mean no limitation like for Arcane
219 if (zi->m_mesh_partitionner_base->nbCellWeight() <= wgt_dim) { // We can use all criteria, so we do
220 cells_weights = zi->m_mesh_partitionner_base->cellsWeightsWithConstraints(wgt_dim);
221 }
222 else { // We need more criteria than available
223 // So we try to balance memory !
224 cells_weights = zi->m_mesh_partitionner_base->cellsSizeWithConstraints();
225 }
226 ArrayView<float> view_weights(cells_weights.size(), obj_weights);
227 view_weights.copy(cells_weights);
228 }
229
230 int index = 0;
231 ENUMERATE_CELL (icell, zi->m_own_cells) {
232 if (!zi->m_mesh_partitionner_base->cellUsedWithConstraints(*icell))
233 continue;
234
235 gids[index] = (*icell).uniqueId().asInt32();
236 lids[index] = zi->m_mesh_partitionner_base->localIdWithConstraints(*icell);
237
238 if (zi->m_ofile) {
239 for (int j = 0; j < wgt_dim; ++j) {
240 float weight = cells_weights[lids[index] * wgt_dim + j];
241 *obj_weights++ = weight;
242 // zi->info() << " Weight uid=" << gids[index] << " w=" << weight;
243 (*zi->m_ofile) << " Weight uid=" << gids[index] << " w=" << weight << '\n';
244 }
245 }
246 index++;
247 }
248 *ierr = ZOLTAN_OK;
249 }
250
251 static void getHgSizeAndFormat(void* data, int* num_lists, int* num_pins, int* format, int* ierr)
252 {
253 ZoltanInfo* zi = (ZoltanInfo*)data;
254
255 zi->info() << " ZoltanInfo::getHgSizeAndFormat() ";
256
257 /*
258 * Supply this query function to Zoltan with Zoltan_Set_HG_Size_CS_Fn().
259 * It tells Zoltan the number of rows or columns to be supplied by
260 * the process, the number of pins (non-zeroes) in the rows or columns,
261 * and whether the pins are provided in compressed row format or
262 * compressed column format.
263 */
264 *format = ZOLTAN_COMPRESSED_EDGE;
265 *num_pins = zi->m_nbPins;
266 *num_lists = zi->m_nbEdges;
267
268 // if (zi->m_model == ZOLTAN_MY_HG) {
269 // *format = ZOLTAN_COMPRESSED_VERTEX;
270 // }
271 zi->info() << " ZoltanInfo::getHgSizeAndFormat() " << " num_list= " << *num_lists << " num_pins= " << *num_pins;
272 *ierr = ZOLTAN_OK;
273 }
274
275 static void getHg(void* data, int num_gid_entries,
276 int nrowcol, int npins, int format,
277 ZOLTAN_ID_PTR z_vtxedge_GID, int* z_vtxedge_ptr, ZOLTAN_ID_PTR z_pin_GID, int* ierr)
278 {
279 ZoltanInfo* zi = (ZoltanInfo*)data;
280
281 zi->info() << " ZoltanInfo::getHg() "
282 << " num_gid_entries= " << num_gid_entries
283 << " nrowcol= " << nrowcol
284 << " npins= " << npins
285 << " format= " << format;
286
287 // 6 neighbors max ? >> Not true for groups !
288 Int64UniqueArray neighbour_cells;
289 neighbour_cells.reserve(6);
290 int indexPin = 0;
291 int indexEdge = 0;
292 int gid = zi->m_edgeGIDStart;
293 z_vtxedge_ptr[0] = 0;
294 ENUMERATE_CELL (i_item, zi->m_own_cells) {
295 const Cell& item = *i_item;
296
297 if (!zi->m_mesh_partitionner_base->cellUsedWithConstraints(item))
298 continue;
299
300 neighbour_cells.resize(0);
301 zi->m_mesh_partitionner_base->getNeighbourCellsUidWithConstraints(item,
302 neighbour_cells);
303
304 if (zi->m_model & ZOLTAN_HG) {
305 for (Integer z = 0; z < neighbour_cells.size(); ++z)
306 z_pin_GID[indexPin++] = CheckedConvert::toInteger(neighbour_cells[z]);
307 z_pin_GID[indexPin++] = item.uniqueId().asInt32(); // Don't forget current cell
308 z_vtxedge_GID[indexEdge] = gid++;
309 // +1 because current cell is it's own neighbor !
310 z_vtxedge_ptr[indexEdge + 1] = z_vtxedge_ptr[indexEdge] + neighbour_cells.size() + 1;
311 indexEdge++;
312 }
313
314 if (zi->m_model & ZOLTAN_GRAPH) {
315 // size 2 edges are face communications
316 // other are neighborhood description
317 for (Integer z = 0; z < neighbour_cells.size(); ++z) {
318 z_pin_GID[indexPin++] = CheckedConvert::toInteger(neighbour_cells[z]);
319 z_pin_GID[indexPin++] = (item.uniqueId().asInt32());
320 z_vtxedge_GID[indexEdge] = gid++;
321 z_vtxedge_ptr[indexEdge + 1] = z_vtxedge_ptr[indexEdge] + 2;
322 indexEdge++;
323 }
324 }
325 }
326 if (zi->m_ofile) {
327 for (int i = 0; i < nrowcol; ++i) {
328 (*zi->m_ofile) << "*topo GID " << z_vtxedge_GID[i]
329 << " : " << z_vtxedge_ptr[i + 1] - z_vtxedge_ptr[i];
330 for (int j = z_vtxedge_ptr[i]; j < z_vtxedge_ptr[i + 1]; ++j) {
331 (*zi->m_ofile) << '\t' << z_pin_GID[j];
332 }
333 (*zi->m_ofile) << '\n';
334 }
335 }
336 *ierr = ZOLTAN_OK;
337 }
338
339 /*
340 For a list of objects, it returns the per-objects sizes (in bytes)
341 of the data buffers needed to pack object data.
342 void (*)(void*, int, int, int, ZOLTAN_ID_TYPE*, ZOLTAN_ID_TYPE*, int*, int*)
343 */
344
345 static void getHgVertexSizes(void* data, int num_gid_entries, int num_lid_entries,
346 int num_ids, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids,
347 int* sizes, int* ierr)
348 {
349 ARCANE_UNUSED(num_gid_entries);
350 ARCANE_UNUSED(num_lid_entries);
351 ARCANE_UNUSED(global_ids);
352
353 ZoltanInfo* zi = (ZoltanInfo*)data;
354
355 SharedArray<float> cellSizes = zi->m_mesh_partitionner_base->cellsSizeWithConstraints();
356
357 for (int i = 0; i < num_ids; ++i) {
358 sizes[i] = Convert::toInteger(cellSizes[local_ids[i]]);
359 }
360
361 *ierr = ZOLTAN_OK;
362 }
363
364 static void getHgEdgeWeightSize(void* data, int* num_edges, int* ierr)
365 {
366 ZoltanInfo* zi = (ZoltanInfo*)data;
367
368 /*
369 * Supply this query function to Zoltan with Zoltan_Set_HG_Size_Edge_Weights_Fn().
370 * It tells Zoltan the number edges for which this process will supply
371 * edge weights. The number of weights per edge was specified with the
372 * parameter EDGE_WEIGHT_DIM. If more than one process provides a weight
373 * for the same edge, the multiple weights are resolved according the
374 * value for the PHG_EDGE_WEIGHT_OPERATION parameter.
375 */
376
377 *num_edges = zi->m_nbEdges;
378 *ierr = ZOLTAN_OK;
379 }
380
381 static void getHgEdgeWeights(void* data, int num_gid_entries,
382 int num_lid_entries, int nedges, int edge_weight_dim,
383 ZOLTAN_ID_PTR edge_GID, ZOLTAN_ID_PTR edge_LID, float* edge_weight, int* ierr)
384 {
385 ARCANE_UNUSED(num_lid_entries);
386
387 ZoltanInfo* zi = (ZoltanInfo*)data;
388 /*
389 * Supply this query function to Zoltan with Zoltan_Set_HG_Edge_Weights_Fn().
390 * It tells Zoltan the weights for some subset of edges.
391 */
392
393 zi->info() << " ZoltanInfo::getHgEdgeWeights() "
394 << " num_gid_entries= " << num_gid_entries
395 << " nedges= " << nedges
396 << " edge_weight_dim= " << edge_weight_dim;
397
398 UniqueArray<float> connectivityWeights;
399 Int64UniqueArray neighbour_cells;
400 neighbour_cells.reserve(6);
401 connectivityWeights.reserve(6);
402 int indexEdge = 0;
403
404 ENUMERATE_CELL (i_item, zi->m_own_cells) {
405 Cell item = *i_item;
406
407 if (!zi->m_mesh_partitionner_base->cellUsedWithConstraints(item))
408 continue;
409
410 connectivityWeights.resize(0);
411 neighbour_cells.resize(0);
412 bool hg_model = (zi->m_model & ZOLTAN_HG);
413 Real he_weight = 0;
414 he_weight =
415 zi->m_mesh_partitionner_base->getNeighbourCellsUidWithConstraints(item,
416 neighbour_cells, &connectivityWeights, hg_model);
417
418 if (zi->m_model & ZOLTAN_HG) {
419 if (!(zi->m_model & ZOLTAN_GRAPH)) {
420 // We have to sum-up de Weight of pins to define HyperEdge's one.
421 for (Integer z = 0; z < neighbour_cells.size(); ++z) {
422 he_weight += connectivityWeights[z];
423 }
424 }
425 edge_GID[indexEdge] = indexEdge + zi->m_edgeGIDStart;
426 edge_LID[indexEdge] = indexEdge;
427 edge_weight[indexEdge++] = (float)he_weight;
428 }
429
430 if (zi->m_model & ZOLTAN_GRAPH) {
431 // size 2 edges are face communications
432 // other are neighborhood description
433 for (Integer z = 0; z < neighbour_cells.size(); ++z) {
434 edge_GID[indexEdge] = indexEdge + zi->m_edgeGIDStart;
435 edge_LID[indexEdge] = indexEdge;
436 Real w = connectivityWeights[z];
437 edge_weight[indexEdge] = (float)w;
438 if (zi->m_mesh_partitionner_base->haveWeakConstraints()) {
439 std::pair<Int64, Int64> items(item.uniqueId(), neighbour_cells[z]);
440 if (zi->m_mesh_partitionner_base->cellUsedWithWeakConstraints(items)) {
441 edge_weight[indexEdge] *= (float)zi->m_edge_weight_multiplier;
442 }
443 }
444 ++indexEdge;
445 }
446 }
447 }
448 *ierr = ZOLTAN_OK;
449 }
450
451 // Return the dimension of a vertex, for geometric methods
452 static int get_num_geometry(void* data, int* ierr)
453 {
454 ARCANE_UNUSED(data);
455 *ierr = ZOLTAN_OK;
456 return 3;
457 }
458
459 // Return the coordinates of my objects (vertices), for geometric methods.
460 static void get_geometry_list(void* data, int sizeGID, int sizeLID,
461 int num_obj,
462 ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids,
463 int num_dim, double* geom_vec, int* ierr)
464 {
465 ARCANE_UNUSED(global_ids);
466 ARCANE_UNUSED(local_ids);
467
468 ZoltanInfo* zi = (ZoltanInfo*)data;
469
470 if ((sizeGID != 1) || (sizeLID != 1) || (num_dim > 3)) {
471 *ierr = ZOLTAN_FATAL;
472 return;
473 }
474
475 VariableNodeReal3& coords(zi->m_mesh_partitionner_base->mesh()->nodesCoordinates());
476
477 int i = 0;
478 ENUMERATE_CELL (icell, zi->m_own_cells) {
479 if (!zi->m_mesh_partitionner_base->cellUsedWithConstraints(*icell))
480 continue;
481
482 // we calculate a barycenter
483
484 Real3 bar;
485
486 for (Integer z = 0, zs = (*icell).nbNode(); z < zs; ++z) {
487 const Node& node = (*icell).node(z);
488 bar += coords[node];
489 }
490 bar /= Convert::toDouble((*icell).nbNode());
491
492 geom_vec[num_dim * i] = bar.x;
493 geom_vec[num_dim * i + 1] = bar.y;
494 geom_vec[num_dim * i + 2] = bar.z;
495 i += 1;
496 }
497
498 String s = platform::getEnvironmentVariable("ZOLTAN_MODEL");
499 if (!s.null() && (s == "MYGEOM")) {
500 for (int i = 0; i < num_obj; ++i) {
501 geom_vec[num_dim * i + num_dim - 1] = 0.0;
502 }
503 }
504
505 *ierr = ZOLTAN_OK;
506 return;
507 }
508};
509
510/*---------------------------------------------------------------------------*/
511/*---------------------------------------------------------------------------*/
512
516class ZoltanMeshPartitioner
518{
519 public:
520
521 explicit ZoltanMeshPartitioner(const ServiceBuildInfo& sbi);
522
523 public:
524
525 virtual void build() {}
526
527 public:
528
529 virtual void partitionMesh(bool initial_partition);
530 virtual void partitionMesh(bool initial_partition, Int32 nb_part);
531
532 virtual void notifyEndPartition();
533
534 private:
535};
536
537/*---------------------------------------------------------------------------*/
538/*---------------------------------------------------------------------------*/
539
540ZoltanMeshPartitioner::
541ZoltanMeshPartitioner(const ServiceBuildInfo& sbi)
542: ArcaneZoltanMeshPartitionerObject(sbi)
543{
544}
545
546/*---------------------------------------------------------------------------*/
547/*---------------------------------------------------------------------------*/
548
550partitionMesh(bool initial_partition)
551{
552 Int32 nb_part = mesh()->parallelMng()->commSize();
553 partitionMesh(initial_partition, nb_part);
554}
555
556/*---------------------------------------------------------------------------*/
557/*---------------------------------------------------------------------------*/
558
560partitionMesh(bool initial_partition, Int32 nb_part)
561{
562 info() << "Load balancing with Zoltan\n";
563
564 float ver;
565
566 int rc = ::Zoltan_Initialize(0, 0, &ver);
567 Zoltan_Memory_Debug(2);
568 if (rc != ZOLTAN_OK)
569 fatal() << "Can not initialize zoltan (r=" << rc << ")";
570 IMesh* mesh = this->mesh();
571 IParallelMng* pm = mesh->parallelMng();
572 Integer nb_rank = pm->commSize();
573
574 if (nb_part < nb_rank)
575 throw ArgumentException(A_FUNCINFO, "partition with nb_part<nb_rank");
576
577 if (nb_rank == 1) {
578 warning() << "Unable to test load balancing on a single sub-domain";
579 return;
580 }
581
582 Integer nb_weight = nbCellWeight();
583
584 struct Zoltan_Struct* zz;
585 int changes;
586 int numGidEntries;
587 int numLidEntries;
588 int numImport;
589 ZOLTAN_ID_PTR importGlobalIds;
590 ZOLTAN_ID_PTR importLocalIds;
591 int* importProcs;
592 int* importToPart;
593 int numExport;
594 ZOLTAN_ID_PTR exportGlobalIds;
595 int* exportProcs;
596
597 /******************************************************************
598 ** Prepare to partition the example hypergraph using the Zoltan
599 ** parallel hypergraph package.
600 **
601 ** Set parameter values, and supply to Zoltan the query functions
602 ** defined in the example library which will provide to the Zoltan
603 ** library the hypergraph data.
604 ******************************************************************/
605 zz = Zoltan_Create(*(MPI_Comm*)getCommunicator());
606
607 Zoltan_Set_Param(zz, "RCB_REUSE", "1"); // Try to do "repartitioning"
608 // option between PHG and RCB (hypergraph and bissection coordinate)
609
610 // Old default values if no environment variable (cf OLD_HG)
611 bool usePHG = true;
612
613 String s = platform::getEnvironmentVariable("ZOLTAN_MODEL");
614
615 if (options()) {
616 if (!s.null() && options()->model.isPresent())
617 fatal() << "Conflicting configuration between ZOLTAN_MODEL environment variable and user data set";
618 usePHG = options()->useHypergraphe();
619 if (!usePHG)
620 s = "RCB";
621 else if (s.null())
622 s = options()->model();
623 }
624 else if (s.null()) {
625 s = "OLDHG";
626 }
627 // => s != null
628
629 if (s == "HYBRID") {
630 if (subDomain()->commonVariables().globalIteration() <= 2)
631 s = "RCB";
632 else
633 s = "DUALHG";
634 }
635
636 String algo = "HYPERGRAPH";
637 if (s == "MYHG" || s == "OLDHG" /* nuance with OLD_HG? */ || s == "DUALHG" || s == "GRAPH") {
638 algo = "HYPERGRAPH";
639 usePHG = true;
640 }
641 else if (s == "RIB" || s == "HSFC") {
642 algo = s;
643 usePHG = false;
644 }
645 else if (s == "RCB" || s == "MYGEOM") {
646 algo = "RCB";
647 usePHG = false;
648 }
649 else {
650 fatal() << "Undefined zoltan model '" << s << "'";
651 }
652
653 int model = ZOLTAN_HG;
654 if (usePHG == false) {
655 model = ZOLTAN_GEOM;
656 }
657
658 if (s == "DUALHG") {
659 model |= ZOLTAN_GRAPH;
660 }
661 else if (s == "GRAPH") {
662 model = ZOLTAN_GRAPH;
663 }
664
665 if (usePHG) {
666 nb_weight = 1; // up to 2 weights for RCB
667 }
668
669 Real edgeWeightMultiplier = 1;
670 Integer repartitionFrequency = 10;
671 Real imbalanceTol = 1.05;
672 Real phgRepartMultiplier = 10;
673 Integer phgOutputLevel = 0;
674 Integer debugLevel = 0;
675
676 if (options()) {
677 edgeWeightMultiplier = options()->edgeWeightMultiplier();
678 repartitionFrequency = options()->repartFrequency();
679 imbalanceTol = options()->imbalanceTol();
680 phgRepartMultiplier = options()->phgRepartMultiplier();
681 phgOutputLevel = options()->phgOutputLevel();
682 debugLevel = options()->debugLevel();
683 }
684
685 /* General parameters */
686 info() << "Zoltan: utilise un repartitionnement " << algo << " (" << s << ").";
687 Zoltan_Set_Param(zz, "LB_METHOD", algo.localstr()); /* partitioning method */
688 Zoltan_Set_Param(zz, "HYPERGRAPH_PACKAGE", "PHG"); /* version of method */
689 //if(!initial_partition)
690 if (mesh->subDomain()->commonVariables().globalIteration() == 1) {
691 Zoltan_Set_Param(zz, "LB_APPROACH", "PARTITION"); /* version of method */
692 info() << "Zoltan: Partition";
693 }
694 else if (mesh->subDomain()->commonVariables().globalIteration() % repartitionFrequency == 0 && repartitionFrequency != -1) {
695 Zoltan_Set_Param(zz, "LB_APPROACH", "REPARTITION"); /* version of method */
696 info() << "Zoltan: Repartition";
697 }
698 else {
699 Zoltan_Set_Param(zz, "LB_APPROACH", "REFINE"); /* version of method */
700 info() << "Zoltan: Refine";
701 }
702
703 String s_imbalaceTol(String::fromNumber(imbalanceTol));
704 Zoltan_Set_Param(zz, "IMBALANCE_TOL", s_imbalaceTol.localstr()); /* version of method */
705
706 String s_nb_part(String::fromNumber(nb_part));
707 Zoltan_Set_Param(zz, "NUM_GLOBAL_PARTS", s_nb_part.localstr());
708
709 Zoltan_Set_Param(zz, "NUM_GID_ENTRIES", "1"); /* global IDs are integers */
710 Zoltan_Set_Param(zz, "NUM_LID_ENTRIES", "1"); /* local IDs are integers */
711 Zoltan_Set_Param(zz, "RETURN_LISTS", "EXPORT"); /* only export lists */
712 String s_nb_weight(String::fromNumber(nb_weight));
713
714 Zoltan_Set_Param(zz, "OBJ_WEIGHT_DIM", s_nb_weight.localstr()); /* 2 vtx weight */
715
716 /* Graph parameters */
717 Zoltan_Set_Param(zz, "ADD_OBJ_WEIGHT", "NONE"); /* Don't calculate extra weights */
718
719 if (usePHG) {
720 /* Parallel hypergraph parameters */
721 Zoltan_Set_Param(zz, "FINAL_OUTPUT", "0"); /* provide stats at the end */
722 Zoltan_Set_Param(zz, "PHG_USE_TIMERS", "1");
723
724 String s_phgOutputLevel(String::fromNumber(phgOutputLevel));
725 Zoltan_Set_Param(zz, "PHG_OUTPUT_LEVEL", s_phgOutputLevel.localstr());
726
727 // Usage in debug
728 Zoltan_Set_Param(zz, "CHECK_HYPERGRAPH", "0"); /* see User's Guide */
729 String s_debugLevel(String::fromNumber(debugLevel));
730 Zoltan_Set_Param(zz, "DEBUG_LEVEL", s_debugLevel.localstr());
731
732 // The default method (ipm) gives theoretically better results
733 // but it seems to sometimes block with Zoltan version 2.0
734
735 Zoltan_Set_Param(zz, "PHG_COARSENING_METHOD", "IPM");
736 //Zoltan_Set_Param(zz, "PHG_COARSEPARTITION_METHOD", "GREEDY");
737 //Zoltan_Set_Param(zz, "PHG_CUT_OBJECTIVE", "HYPEREDGES");
738
739 Zoltan_Set_Param(zz, "EDGE_WEIGHT_DIM", "1");
740 Zoltan_Set_Param(zz, "PHG_EDGE_WEIGHT_OPERATION", "add");
741
742 if (!initial_partition) {
743 // Number of iterations between 2 load balances. Used to scale the cost
744 // of data migration.
745 String s_phgRepartMultiplier(String::fromNumber(phgRepartMultiplier));
746 Zoltan_Set_Param(zz, "PHG_REPART_MULTIPLIER", s_phgRepartMultiplier.localstr());
747 }
748 }
749 else {
750 Zoltan_Set_Param(zz, "KEEP_CUTS", "0");
751 Zoltan_Set_Param(zz, "RCB_OUTPUT_LEVEL", "0");
752 Zoltan_Set_Param(zz, "RCB_RECTILINEAR_BLOCKS", "0");
753 // The 2 following options are for multi-weights partitioning
754 // Object weights are not comparable
755 Zoltan_Set_Param(zz, "OBJ_WEIGHTS_COMPARABLE", "0");
756 // 1 : minimize total (ie Sum over all phases)
757 // 2 : between
758 // 3 : Minimize worst case for each phase
759 Zoltan_Set_Param(zz, "RCB_MULTICRITERIA_NORM", "3");
760 Zoltan_Set_Param(zz, "RCB_MAX_ASPECT_RATIO", "10");
761 }
762
763 bool dump_infos = false;
764 if (platform::getEnvironmentVariable("ARCANE_DEBUG_PARTITION") == "TRUE")
765 dump_infos = true;
766 ofstream ofile;
767 if (dump_infos) {
768 StringBuilder fname;
769 Integer iteration = mesh->subDomain()->commonVariables().globalIteration();
770 fname = "weigth-";
771 fname += pm->commRank();
772 fname += "-";
773 fname += iteration;
774 String f(fname);
775 ofile.open(f.localstr());
776 }
777
778 /* Application defined query functions (defined in exphg.c) */
779
780 // initialisation pour la gestion des contraintes
781 initConstraints();
782
783 ostream* zofile = 0;
784 if (dump_infos)
785 zofile = &ofile;
786
787 ScopedPtrT<ZoltanInfo> zoltan_info(new ZoltanInfo(this, zofile, model, edgeWeightMultiplier));
788 Zoltan_Set_Num_Obj_Fn(zz, &ZoltanInfo::getHgNumVertices, zoltan_info.get());
789 Zoltan_Set_Obj_List_Fn(zz, &ZoltanInfo::getHgVerticesAndWeights, zoltan_info.get());
790 if (!initial_partition) {
791 Zoltan_Set_Obj_Size_Multi_Fn(zz, &ZoltanInfo::getHgVertexSizes,
792 zoltan_info.get());
793 }
794 if (usePHG) {
795 info() << "Setting up HG Callbacks";
796 Zoltan_Set_HG_Size_CS_Fn(zz, &ZoltanInfo::getHgSizeAndFormat, zoltan_info.get());
797 Zoltan_Set_HG_CS_Fn(zz, &ZoltanInfo::getHg, zoltan_info.get());
798 Zoltan_Set_HG_Size_Edge_Weights_Fn(zz, &ZoltanInfo::getHgEdgeWeightSize, zoltan_info.get());
799 Zoltan_Set_HG_Edge_Weights_Fn(zz, &ZoltanInfo::getHgEdgeWeights, zoltan_info.get());
800 }
801 else {
802 info() << "Setting up Geom Callbacks";
803 Zoltan_Set_Num_Geom_Fn(zz, ZoltanInfo::get_num_geometry, zoltan_info.get());
804 Zoltan_Set_Geom_Multi_Fn(zz, ZoltanInfo::get_geometry_list, zoltan_info.get());
805 }
806
807 /* Parallel partitioning occurs now */
808
809 int* export_partitions = 0;
810 ZOLTAN_ID_PTR export_local_ids = 0;
811
812 info() << "Doing partition";
813 rc = Zoltan_LB_Partition(zz, &changes, &numGidEntries, &numLidEntries,
814 &numImport, &importGlobalIds, &importLocalIds,
815 &importProcs, &importToPart,
816 &numExport, &exportGlobalIds, &export_local_ids,
817 &exportProcs, &export_partitions);
818 pm->barrier();
819
820 VariableItemInt32& cells_new_owner = mesh->toPrimaryMesh()->itemsNewOwner(IK_Cell);
821 ENUMERATE_ITEM (icell, mesh->ownCells()) {
822 cells_new_owner[icell] = (*icell).owner();
823 }
824
825 invertArrayLid2LidCompacted();
826
827 {
828 CellInfoListView items_internal(m_cell_family);
829 Integer nb_export = numExport;
830 info() << "numExport = " << numExport;
831 for (Integer i = 0; i < nb_export; ++i) {
832 Item item = items_internal[localIdWithConstraints(export_local_ids[i])];
833 // change for the cell or the whole group if necessary
834 changeCellOwner(item, cells_new_owner, export_partitions[i]);
835 if (dump_infos) // these infos do not take groups/constraints into account
836 ofile << "EXPORT: uid=" << ItemPrinter(item) << " OLD=" << item.owner()
837 << " NEW=" << cells_new_owner[item] << " PROC=" << exportProcs[i] << '\n';
838 }
839 // pinfo() << "Proc " << my_rank << " number of cells changing domain: "
840 // << nb_export << " changes=" << changes
841 // << " ZoltanMem=" << Zoltan_Memory_Usage(ZOLTAN_MEM_STAT_MAXIMUM);
842 }
843
844 // Free the memory of zoltan_info
845 zoltan_info = 0;
846 if (dump_infos)
847 ofile.close();
848
849 if (numExport < 0) {
850 sleep(3);
851 MPI_Abort(MPI_COMM_WORLD, EXIT_FAILURE);
852 }
853
854 Zoltan_LB_Free_Part(&exportGlobalIds, &export_local_ids,
855 &exportProcs, &export_partitions);
856
857 // Zoltan_Destroy(&zz);
858
859 // release of temporary arrays
860 freeConstraints();
861
862 cells_new_owner.synchronize();
864}
865
866/*---------------------------------------------------------------------------*/
867/*---------------------------------------------------------------------------*/
868
873
874/*---------------------------------------------------------------------------*/
875/*---------------------------------------------------------------------------*/
876
878 ServiceProperty("Zoltan", ST_SubDomain),
881ARCANE_REGISTER_SERVICE_ZOLTANMESHPARTITIONER(Zoltan, ZoltanMeshPartitioner);
882
883#if ARCANE_DEFAULT_PARTITIONER == ZOLTAN_DEFAULT_PARTITIONER
885 ServiceProperty("DefaultPartitioner", ST_SubDomain),
888ARCANE_REGISTER_SERVICE_ZOLTANMESHPARTITIONER(DefaultPartitioner, ZoltanMeshPartitioner);
889#endif
890
891/*---------------------------------------------------------------------------*/
892/*---------------------------------------------------------------------------*/
893
894} // End namespace Arcane
895
896/*---------------------------------------------------------------------------*/
897/*---------------------------------------------------------------------------*/
Types and macros for iterating over mesh entities.
#define ENUMERATE_FACE(name, group)
Generic enumerator for a face group.
#define ENUMERATE_CELL(name, group)
Generic enumerator for a cell group.
#define ENUMERATE_ITEM(name, group)
Generic enumerator for a node group.
#define ARCANE_SERVICE_INTERFACE(ainterface)
Macro to declare an interface when registering a service.
Integer size() const
Number of elements in the vector.
CaseOptionsZoltanMeshPartitioner * options() const
Options du jeu de données du service.
ArcaneZoltanMeshPartitionerObject(const Arcane::ServiceBuildInfo &sbi)
Constructeur.
Modifiable view of an array of type T.
void copy(const U &copy_array)
Copies the array copy_array into the instance.
void fill(ConstReferenceType value)
Fills the array with the value value.
void resize(Int64 s)
Changes the number of elements in the array to s.
void reserve(Int64 new_capacity)
Reserves memory for new_capacity elements.
Cell of a mesh.
Definition Item.h:1300
Face of a cell.
Definition Item.h:1032
Cell frontCell() const
Cell in front of the face (null cell if none).
Definition Item.h:1780
Cell backCell() const
Cell behind the face (null cell if none).
Definition Item.h:1774
virtual ITraceMng * traceMng()=0
Associated message manager.
Interface of a mesh partitioner.
Interface of a mesh partitioner.
virtual IMesh * mesh() const =0
Mesh associated with the partitioner.
virtual VariableNodeReal3 & nodesCoordinates()=0
Node coordinates.
Interface of the parallelism manager for a subdomain.
virtual Int32 commRank() const =0
Rank of this instance in the communicator.
virtual void scan(eReduceType rt, ArrayView< char > v)=0
Applies a prefix-sum algorithm on the values of v using the rt operation.
virtual Int32 commSize() const =0
Number of instances in the communicator.
virtual void barrier()=0
Performs a barrier.
ItemUniqueId uniqueId() const
Unique identifier across all domains.
Definition Item.h:239
constexpr bool null() const
true if the entity is null (i.e. not connected to the mesh)
Definition Item.h:230
Base class for a load balancing service.
IMesh * mesh() const override
Mesh associated with the partitioner.
virtual void changeOwnersFromCells()
Positions the new owners of nodes, edges and faces based on the cells.
Scalar variable on a mesh entity type.
Node of a mesh.
Definition Item.h:598
Class managing a 3-dimensional real vector.
Definition Real3.h:132
Structure containing the information to create a service.
Service creation properties.
1D vector of data with reference semantics.
bool null() const
Returns true if the string is null.
Definition String.cc:306
TraceAccessor(ITraceMng *m)
Constructs an accessor via the trace manager m.
TraceMessage fatal() const
Flow for a fatal error message.
TraceMessage info() const
Flow for an information message.
TraceMessage warning() const
Flow for a warning message.
1D data vector with value semantics (STL style).
Parameters necessary for building a variable.
Mesh partitioner using the Zoltan library.
virtual void notifyEndPartition()
Notification when a re-partitioning finishes (after entity exchange).
virtual void build()
Build-level construction of the service.
virtual void partitionMesh(bool initial_partition)
ItemGroupT< Cell > CellGroup
Group of cells.
Definition ItemTypes.h:184
#define ARCANE_REGISTER_SERVICE(aclass, a_service_property,...)
Macro for registering a service.
MeshVariableScalarRefT< Node, Real3 > VariableNodeReal3
Coordinate type quantity at node.
ItemVariableScalarRefT< Int32 > VariableItemInt32
32-bit integer type quantity
Integer toInteger(Real r)
Converts a Real to Integer.
double toDouble(Real r)
Converts a Real to double.
void sleep(Integer nb_second)
Puts the process to sleep for nb_second seconds.
String getEnvironmentVariable(const String &name)
Environment variable named name.
-- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature --
UniqueArray< Int64 > Int64UniqueArray
Dynamic 1D array of 64-bit integers.
Definition UtilsTypes.h:339
Int32 Integer
Type representing an integer.
@ ST_SubDomain
The service is used at the subdomain level.
UniqueArray< Int32 > Int32UniqueArray
Dynamic 1D array of 32-bit integers.
Definition UtilsTypes.h:341
@ IK_Cell
Cell mesh entity.
double Real
Type representing a real number.
std::int32_t Int32
Signed integer type of 32 bits.